Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4EmDNAPhysics.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// $Id$
27
28#include "G4EmDNAPhysics.hh"
29
30#include "G4SystemOfUnits.hh"
31
33
34// *** Processes and models for Geant4-DNA
35
36#include "G4DNAElastic.hh"
39
40#include "G4DNAExcitation.hh"
41#include "G4DNAAttachment.hh"
42#include "G4DNAVibExcitation.hh"
43#include "G4DNAIonisation.hh"
46
47// particles
48
49#include "G4Electron.hh"
50#include "G4Proton.hh"
51#include "G4GenericIon.hh"
52
53// Warning : the following is needed in order to use EM Physics builders
54// e+
55#include "G4Positron.hh"
57#include "G4UrbanMscModel95.hh"
58#include "G4eIonisation.hh"
59#include "G4eBremsstrahlung.hh"
61// gamma
62#include "G4Gamma.hh"
67#include "G4GammaConversion.hh"
71// end of warning
72
73#include "G4LossTableManager.hh"
76#include "G4BuilderType.hh"
77
78// factory
80//
82
83
84//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
85
87 : G4VPhysicsConstructor("G4EmDNAPhysics"), verbose(ver)
88{
90}
91
92//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
93
95 : G4VPhysicsConstructor("G4EmDNAPhysics"), verbose(ver)
96{
98}
99
100//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
101
103{}
104
105//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
106
108{
109// bosons
111
112// leptons
115
116// baryons
118
120
121 G4DNAGenericIonsManager * genericIonsManager;
122 genericIonsManager=G4DNAGenericIonsManager::Instance();
123 genericIonsManager->GetIon("alpha++");
124 genericIonsManager->GetIon("alpha+");
125 genericIonsManager->GetIon("helium");
126 genericIonsManager->GetIon("hydrogen");
127 genericIonsManager->GetIon("carbon");
128 genericIonsManager->GetIon("nitrogen");
129 genericIonsManager->GetIon("oxygen");
130 genericIonsManager->GetIon("iron");
131
132}
133
134//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
135
137{
139
141 while( (*theParticleIterator)() )
142 {
144 G4String particleName = particle->GetParticleName();
145
146 if (particleName == "e-") {
147
148 // *** Elastic scattering (two alternative models available) ***
149
150 G4DNAElastic* theDNAElasticProcess = new G4DNAElastic("e-_G4DNAElastic");
151 theDNAElasticProcess->SetEmModel(new G4DNAChampionElasticModel());
152
153 // or alternative model
154 //theDNAElasticProcess->SetEmModel(new G4DNAScreenedRutherfordElasticModel());
155
156 ph->RegisterProcess(theDNAElasticProcess, particle);
157
158 // *** Excitation ***
159 ph->RegisterProcess(new G4DNAExcitation("e-_G4DNAExcitation"), particle);
160
161 // *** Ionisation ***
162 ph->RegisterProcess(new G4DNAIonisation("e-_G4DNAIonisation"), particle);
163
164 // *** Vibrational excitation ***
165 ph->RegisterProcess(new G4DNAVibExcitation("e-_G4DNAVibExcitation"), particle);
166
167 // *** Attachment ***
168 ph->RegisterProcess(new G4DNAAttachment("e-_G4DNAAttachment"), particle);
169
170 } else if ( particleName == "proton" ) {
171 ph->RegisterProcess(new G4DNAExcitation("proton_G4DNAExcitation"), particle);
172 ph->RegisterProcess(new G4DNAIonisation("proton_G4DNAIonisation"), particle);
173 ph->RegisterProcess(new G4DNAChargeDecrease("proton_G4DNAChargeDecrease"), particle);
174
175 } else if ( particleName == "hydrogen" ) {
176 ph->RegisterProcess(new G4DNAExcitation("hydrogen_G4DNAExcitation"), particle);
177 ph->RegisterProcess(new G4DNAIonisation("hydrogen_G4DNAIonisation"), particle);
178 ph->RegisterProcess(new G4DNAChargeIncrease("hydrogen_G4DNAChargeIncrease"), particle);
179
180 } else if ( particleName == "alpha" ) {
181 ph->RegisterProcess(new G4DNAExcitation("alpha_G4DNAExcitation"), particle);
182 ph->RegisterProcess(new G4DNAIonisation("alpha_G4DNAIonisation"), particle);
183 ph->RegisterProcess(new G4DNAChargeDecrease("alpha_G4DNAChargeDecrease"), particle);
184
185 } else if ( particleName == "alpha+" ) {
186 ph->RegisterProcess(new G4DNAExcitation("alpha+_G4DNAExcitation"), particle);
187 ph->RegisterProcess(new G4DNAIonisation("alpha+_G4DNAIonisation"), particle);
188 ph->RegisterProcess(new G4DNAChargeDecrease("alpha+_G4DNAChargeDecrease"), particle);
189 ph->RegisterProcess(new G4DNAChargeIncrease("alpha+_G4DNAChargeIncrease"), particle);
190
191 } else if ( particleName == "helium" ) {
192 ph->RegisterProcess(new G4DNAExcitation("helium_G4DNAExcitation"), particle);
193 ph->RegisterProcess(new G4DNAIonisation("helium_G4DNAIonisation"), particle);
194 ph->RegisterProcess(new G4DNAChargeIncrease("helium_G4DNAChargeIncrease"), particle);
195
196 // Extension to HZE proposed by Z. Francis
197
198 } else if ( particleName == "carbon" ) {
199 ph->RegisterProcess(new G4DNAIonisation("carbon_G4DNAIonisation"), particle);
200
201 } else if ( particleName == "nitrogen" ) {
202 ph->RegisterProcess(new G4DNAIonisation("nitrogen_G4DNAIonisation"), particle);
203
204 } else if ( particleName == "oxygen" ) {
205 ph->RegisterProcess(new G4DNAIonisation("oxygen_G4DNAIonisation"), particle);
206
207 } else if ( particleName == "iron" ) {
208 ph->RegisterProcess(new G4DNAIonisation("iron_G4DNAIonisation"), particle);
209
210 }
211
212 // Warning : the following particles and processes are needed by EM Physics builders
213 // They are taken from the default Livermore Physics list
214 // These particles are currently not handled by Geant4-DNA
215
216 // e+
217
218 else if (particleName == "e+") {
219
220 // Identical to G4EmStandardPhysics_option3
221
223 msc->AddEmModel(0, new G4UrbanMscModel95());
225 G4eIonisation* eIoni = new G4eIonisation();
226 eIoni->SetStepFunction(0.2, 100*um);
227
228 ph->RegisterProcess(msc, particle);
229 ph->RegisterProcess(eIoni, particle);
230 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
231 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
232
233 } else if (particleName == "gamma") {
234
235 G4double LivermoreHighEnergyLimit = GeV;
236
237 G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
238 G4LivermorePhotoElectricModel* theLivermorePhotoElectricModel =
240 theLivermorePhotoElectricModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
241 thePhotoElectricEffect->AddEmModel(0, theLivermorePhotoElectricModel);
242 ph->RegisterProcess(thePhotoElectricEffect, particle);
243
244 G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
245 G4LivermoreComptonModel* theLivermoreComptonModel =
247 theLivermoreComptonModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
248 theComptonScattering->AddEmModel(0, theLivermoreComptonModel);
249 ph->RegisterProcess(theComptonScattering, particle);
250
251 G4GammaConversion* theGammaConversion = new G4GammaConversion();
252 G4LivermoreGammaConversionModel* theLivermoreGammaConversionModel =
254 theLivermoreGammaConversionModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
255 theGammaConversion->AddEmModel(0, theLivermoreGammaConversionModel);
256 ph->RegisterProcess(theGammaConversion, particle);
257
258 G4RayleighScattering* theRayleigh = new G4RayleighScattering();
259 G4LivermoreRayleighModel* theRayleighModel = new G4LivermoreRayleighModel();
260 theRayleighModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
261 theRayleigh->AddEmModel(0, theRayleighModel);
262 ph->RegisterProcess(theRayleigh, particle);
263 }
264
265 // Warning : end of particles and processes are needed by EM Physics builders
266
267 }
268
269 // Deexcitation
270 //
273 de->SetFluo(true);
274}
275
276//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ bElectromagnetic
@ fUseDistanceToBoundary
#define G4_DECLARE_PHYSCONSTR_FACTORY(physics_constructor)
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
static G4DNAGenericIonsManager * Instance(void)
G4ParticleDefinition * GetIon(const G4String &name)
static G4Electron * Electron()
Definition: G4Electron.cc:94
virtual void ConstructProcess()
virtual void ConstructParticle()
virtual ~G4EmDNAPhysics()
G4EmDNAPhysics(G4int ver=1)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4GenericIon * GenericIonDefinition()
Definition: G4GenericIon.cc:87
void SetAtomDeexcitation(G4VAtomDeexcitation *)
static G4LossTableManager * Instance()
const G4String & GetParticleName() const
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4PhysicsListHelper * GetPhysicsListHelper()
static G4Positron * Positron()
Definition: G4Positron.cc:94
static G4Proton * Proton()
Definition: G4Proton.cc:93
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=0)
void SetEmModel(G4VEmModel *, G4int index=1)
void SetStepFunction(G4double v1, G4double v2)
void AddEmModel(G4int order, G4VEmModel *, const G4Region *region=0)
void SetStepLimitType(G4MscStepLimitType val)
G4ParticleTable::G4PTblDicIterator * theParticleIterator