source: trunk/examples/advanced/human_phantom/src/G4HumanPhantomPhysicsList.cc@ 810

Last change on this file since 810 was 807, checked in by garnier, 17 years ago

update

File size: 5.8 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// Authors: S. Guatelli and M. G. Pia, INFN Genova, Italy
27//
28// Based on code developed by the undergraduate student G. Guerrieri
29// Note: this is a preliminary beta-version of the code; an improved
30// version will be distributed in the next Geant4 public release, compliant
31// with the design in a forthcoming publication, and subject to a
32// design and code review.
33
34#include "G4HumanPhantomPhysicsList.hh"
35
36#include "G4ParticleDefinition.hh"
37#include "G4ProductionCutsTable.hh"
38#include "G4ProcessManager.hh"
39#include "G4ParticleTypes.hh"
40#include "G4UnitsTable.hh"
41#include "G4ios.hh"
42
43#include "G4MultipleScattering.hh"
44// gamma
45#include "G4LowEnergyRayleigh.hh"
46#include "G4LowEnergyPhotoElectric.hh"
47#include "G4LowEnergyCompton.hh"
48#include "G4LowEnergyGammaConversion.hh"
49// e-
50#include "G4LowEnergyIonisation.hh"
51#include "G4LowEnergyBremsstrahlung.hh"
52// e+
53#include "G4eIonisation.hh"
54#include "G4eBremsstrahlung.hh"
55#include "G4eplusAnnihilation.hh"
56
57G4HumanPhantomPhysicsList::G4HumanPhantomPhysicsList(): G4VUserPhysicsList()
58{
59 SetVerboseLevel(1);
60}
61
62G4HumanPhantomPhysicsList::~G4HumanPhantomPhysicsList()
63{
64}
65
66void G4HumanPhantomPhysicsList::ConstructParticle()
67{
68 // In this method, static member functions should be called
69 // for all particles which you want to use.
70 // This ensures that objects of these particle types will be
71 // created in the program.
72 G4Geantino::GeantinoDefinition();
73
74 ConstructBosons();
75 ConstructLeptons();
76}
77
78void G4HumanPhantomPhysicsList::ConstructBosons()
79{
80 // photons
81 G4Gamma::GammaDefinition();
82}
83
84void G4HumanPhantomPhysicsList::ConstructLeptons()
85{
86 // leptons
87 G4Electron::ElectronDefinition();
88 G4Positron::PositronDefinition();
89}
90
91void G4HumanPhantomPhysicsList::ConstructProcess()
92{
93 AddTransportation();
94 ConstructEM();
95}
96
97void G4HumanPhantomPhysicsList::ConstructEM()
98{
99 theParticleIterator->reset();
100
101 while( (*theParticleIterator)() ){
102
103 G4ParticleDefinition* particle = theParticleIterator->value();
104 G4ProcessManager* pmanager = particle->GetProcessManager();
105 G4String particleName = particle->GetParticleName();
106
107 // Processes
108
109 if (particleName == "gamma") {
110 // Photon
111 pmanager->AddDiscreteProcess(new G4LowEnergyRayleigh);
112 pmanager->AddDiscreteProcess(new G4LowEnergyPhotoElectric);
113 pmanager->AddDiscreteProcess(new G4LowEnergyCompton);
114 pmanager->AddDiscreteProcess(new G4LowEnergyGammaConversion);
115
116 } else if (particleName == "e-") {
117 // Electron
118 G4LowEnergyIonisation* loweIon = new G4LowEnergyIonisation("LowEnergyIoni");
119
120 G4LowEnergyBremsstrahlung* loweBrem = new G4LowEnergyBremsstrahlung("LowEnBrem");
121 // Select the Bremsstrahlung angular distribution model (Tsai/2BN/2BS)
122 loweBrem->SetAngularGenerator("tsai");
123
124 pmanager->AddProcess(new G4MultipleScattering, -1, 1,1);
125 pmanager->AddProcess(loweIon, -1, 2,2);
126 pmanager->AddProcess(loweBrem, -1,-1,3);
127
128 } else if (particleName == "e+") {
129 // Positron
130 pmanager->AddProcess(new G4MultipleScattering, -1, 1,1);
131 pmanager->AddProcess(new G4eIonisation, -1, 2,2);
132 pmanager->AddProcess(new G4eBremsstrahlung, -1,-1,3);
133 pmanager->AddProcess(new G4eplusAnnihilation, 0,-1,4);
134
135 }
136 }
137}
138
139void G4HumanPhantomPhysicsList::SetCuts()
140{
141 // The production threshold is fixed to 0.1 mm for all the particles
142 // Secondary particles with a range bigger than 0.1 mm
143 // are generated; otherwise their energy is considered deposited locally
144
145 defaultCutValue = 0.1 * mm;
146
147 const G4double cutForGamma = defaultCutValue;
148 const G4double cutForElectron = defaultCutValue;
149 const G4double cutForPositron = defaultCutValue;
150
151 SetCutValue(cutForGamma, "gamma");
152 SetCutValue(cutForElectron, "e-");
153 SetCutValue(cutForPositron, "e+");
154
155 // Set the secondary production cut lower than 990. eV
156 // Very important for high precision of lowenergy processes at low energies
157
158 G4double lowLimit = 250. * eV;
159 G4double highLimit = 100. * GeV;
160 G4ProductionCutsTable::GetProductionCutsTable()->SetEnergyRange(lowLimit, highLimit);
161
162 if (verboseLevel>0) DumpCutValuesTable();
163}
Note: See TracBrowser for help on using the repository browser.