SToGS
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros
SToGS_LowEnergyEMPhysicsList.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 
29 
30 #include "G4ProcessManager.hh"
31 #include "G4ProcessVector.hh"
32 #include "G4LossTableManager.hh"
33 #include "G4ParticleDefinition.hh"
34 #include "G4ParticleWithCuts.hh"
35 #include "G4ParticleTypes.hh"
36 #include "G4ParticleTable.hh"
37 #include "G4ios.hh"
38 
39 // gamma
40 #include "G4PhotoElectricEffect.hh"
41 #include "G4LivermorePhotoElectricModel.hh"
42 
43 #include "G4ComptonScattering.hh"
44 #include "G4LivermoreComptonModel.hh"
45 
46 #include "G4GammaConversion.hh"
47 #include "G4LivermoreGammaConversionModel.hh"
48 
49 #include "G4RayleighScattering.hh"
50 #include "G4LivermoreRayleighModel.hh"
51 
52 // e-
53 #include "G4eMultipleScattering.hh"
54 #include "G4UniversalFluctuation.hh"
55 
56 #include "G4eIonisation.hh"
57 #include "G4LivermoreIonisationModel.hh"
58 
59 #include "G4eBremsstrahlung.hh"
60 #include "G4LivermoreBremsstrahlungModel.hh"
61 
62 #include "G4eIonisation.hh"
63 #include "G4eBremsstrahlung.hh"
64 #include "G4eplusAnnihilation.hh"
65 
66 
67 ParisLowEnergyEMPhysicsList::ParisLowEnergyEMPhysicsList(const G4String& name): G4VPhysicsConstructor(name)
68 {
69 }
70 
71 ParisLowEnergyEMPhysicsList::~ParisLowEnergyEMPhysicsList()
72 {
73 }
74 
75 void ParisLowEnergyEMPhysicsList::ConstructParticle()
76 {
77  // gamma
78  G4Gamma::GammaDefinition();
79 
80  // electron
81  G4Electron::ElectronDefinition();
82  G4Positron::PositronDefinition();
83 
84  // pseudo-particles
85  G4Geantino::GeantinoDefinition();
86 }
87 
88 
89 void ParisLowEnergyEMPhysicsList::ConstructProcess()
90 {
91  theParticleIterator->reset();
92  //
93  while( (*theParticleIterator)() ){
94  G4ParticleDefinition* particle = theParticleIterator->value();
95  G4ProcessManager* pmanager = particle->GetProcessManager();
96  G4String particleName = particle->GetParticleName();
97 
98  if (particleName == "gamma") {
99 
100  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
101  thePhotoElectricEffect->SetModel(new G4LivermorePhotoElectricModel());
102 
103  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
104  theComptonScattering->SetModel(new G4LivermoreComptonModel());
105 
106  G4GammaConversion* theGammaConversion = new G4GammaConversion();
107  theGammaConversion->SetModel(new G4LivermoreGammaConversionModel());
108 
109  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
110  theRayleigh->SetModel(new G4LivermoreRayleighModel());
111 
112  pmanager->AddDiscreteProcess(thePhotoElectricEffect);
113  pmanager->AddDiscreteProcess(theComptonScattering);
114  pmanager->AddDiscreteProcess(theGammaConversion);
115  pmanager->AddDiscreteProcess(theRayleigh);
116 
117  } else if (particleName == "e-") {
118 
119  G4eMultipleScattering* msc = new G4eMultipleScattering();
120 
121  // Ionisation
122  G4eIonisation* eIoni = new G4eIonisation();
123  eIoni->SetEmModel(new G4LivermoreIonisationModel());
124  eIoni->SetFluctModel(new G4UniversalFluctuation() );
125 
126  // Bremsstrahlung
127  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
128  eBrem->SetEmModel(new G4LivermoreBremsstrahlungModel());
129 
130  pmanager->AddProcess(msc ,-1, 1,1);
131  pmanager->AddProcess(eIoni ,-1, 2,2);
132  pmanager->AddProcess(eBrem ,-1,-1,3);
133 
134  } else if (particleName == "e+") {
135 
136  G4eMultipleScattering* msc = new G4eMultipleScattering();
137 
138  // Ionisation
139  G4eIonisation* eIoni = new G4eIonisation();
140  eIoni->SetEmModel(new G4LivermoreIonisationModel());
141  eIoni->SetFluctModel(new G4UniversalFluctuation());
142 
143  // Bremsstrahlung
144  G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
145  eBrem->SetEmModel(new G4LivermoreBremsstrahlungModel());
146 
147  //Annihilation
148  G4eplusAnnihilation* eAnni = new G4eplusAnnihilation();
149 
150  pmanager->AddProcess(msc ,-1, 1,1);
151  pmanager->AddProcess(eIoni ,-1, 2,2);
152  pmanager->AddProcess(eBrem ,-1,-1,3);
153  pmanager->AddProcess(eAnni , 0,-1,4);
154 
155  }
156  }
157 }
158 
159 /*
160 void ParisLowEnergyEMPhysicsList::ConstructEM()
161 {
162 // G4VProcess *aprocess;
163 
164  theParticleIterator->reset();
165  while( (*theParticleIterator)() ){
166  G4ParticleDefinition* particle = theParticleIterator->value();
167  G4ProcessManager* pmanager = particle->GetProcessManager();
168  G4String particleName = particle->GetParticleName();
169 
170  if (particleName == "gamma") {
171  pmanager->AddDiscreteProcess(new G4LowEnergyPhotoElectric);
172  pmanager->AddDiscreteProcess(new G4LowEnergyCompton);
173  pmanager->AddDiscreteProcess(new G4LowEnergyGammaConversion);
174  pmanager->AddDiscreteProcess(new G4LowEnergyRayleigh);
175 
176  } else if (particleName == "e-") {
177  pmanager->AddProcess(new G4eMultipleScattering, -1, 1,1);
178  pmanager->AddProcess(new G4LowEnergyIonisation, -1, 2,2);
179  pmanager->AddProcess(new G4LowEnergyBremsstrahlung, -1,-1,3);
180 
181  } else if (particleName == "e+") {
182  pmanager->AddProcess(new G4eMultipleScattering, -1, 1,1);
183  pmanager->AddProcess(new G4eIonisation, -1, 2,2);
184  pmanager->AddProcess(new G4eBremsstrahlung, -1,-1,3);
185  pmanager->AddProcess(new G4eplusAnnihilation, 0,-1,4);
186 
187  }
188  }
189 }
190 */
191 
192 void ParisLowEnergyEMPhysicsList::SetCuts()
193 {
194  /*
195  if ( verboseLevel > 0 ) {
196  G4cout << " ParisLowEnergyEMPhysicsList::SetCuts:";
197  G4cout << " CutLength : " << G4BestUnit(defaultCutValue,"Length") << G4endl;
198  }
199 
200  // set cut values for gamma at first and for e- second and next for e+,
201  // because some processes for e+/e- need cut values for gamma
202  //
203  SetCutValue(defaultCutValue, "gamma");
204  SetCutValue(defaultCutValue, "e-");
205  SetCutValue(defaultCutValue, "e+");
206 
207  if (verboseLevel>0) DumpCutValuesTable();
208  */
209 }
210 
211 
212 
213