Geant4  10.01
G4EmDNAPhysics_option1.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: G4EmDNAPhysics_option1.cc 73244 2013-08-22 13:12:17Z gcosmo $
27 
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"
44 #include "G4DNAChargeDecrease.hh"
45 #include "G4DNAChargeIncrease.hh"
46 
47 #include "G4hMultipleScattering.hh"
48 #include "G4LowEWentzelVIModel.hh"
49 
50 // particles
51 
52 #include "G4Electron.hh"
53 #include "G4Proton.hh"
54 #include "G4GenericIon.hh"
55 
56 // Warning : the following is needed in order to use EM Physics builders
57 // e+
58 #include "G4Positron.hh"
59 #include "G4eMultipleScattering.hh"
60 #include "G4eIonisation.hh"
61 #include "G4eBremsstrahlung.hh"
62 #include "G4eplusAnnihilation.hh"
63 // gamma
64 #include "G4Gamma.hh"
65 #include "G4PhotoElectricEffect.hh"
67 #include "G4ComptonScattering.hh"
69 #include "G4GammaConversion.hh"
71 #include "G4RayleighScattering.hh"
73 // end of warning
74 
75 #include "G4LossTableManager.hh"
76 #include "G4UAtomicDeexcitation.hh"
77 #include "G4PhysicsListHelper.hh"
78 #include "G4BuilderType.hh"
79 
80 // factory
82 //
84 
85 
86 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
87 
89  : G4VPhysicsConstructor("G4EmDNAPhysics_option1"), verbose(ver)
90 {
92 }
93 
94 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
95 
97  : G4VPhysicsConstructor("G4EmDNAPhysics_option1"), verbose(ver)
98 {
100 }
101 
102 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
103 
105 {}
106 
107 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
108 
110 {
111 // bosons
112  G4Gamma::Gamma();
113 
114 // leptons
117 
118 // baryons
120 
122 
123  G4DNAGenericIonsManager * genericIonsManager;
124  genericIonsManager=G4DNAGenericIonsManager::Instance();
125  genericIonsManager->GetIon("alpha++");
126  genericIonsManager->GetIon("alpha+");
127  genericIonsManager->GetIon("helium");
128  genericIonsManager->GetIon("hydrogen");
129  genericIonsManager->GetIon("carbon");
130  genericIonsManager->GetIon("nitrogen");
131  genericIonsManager->GetIon("oxygen");
132  genericIonsManager->GetIon("iron");
133 
134 }
135 
136 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
137 
139 {
140  if(verbose > 1) {
141  G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
142  }
144 
145  aParticleIterator->reset();
146  while( (*aParticleIterator)() )
147  {
148  G4ParticleDefinition* particle = aParticleIterator->value();
149  G4String particleName = particle->GetParticleName();
150 
151  if (particleName == "e-") {
152 
153  // *** Elastic scattering (two alternative models available) ***
154 
155  //G4DNAElastic* theDNAElasticProcess = new G4DNAElastic("e-_G4DNAElastic");
156  //theDNAElasticProcess->SetEmModel(new G4DNAChampionElasticModel());
157 
158  // or alternative model
159  //theDNAElasticProcess->SetEmModel(new G4DNAScreenedRutherfordElasticModel());
160 
161  //ph->RegisterProcess(theDNAElasticProcess, particle);
162 
164  msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
165  ph->RegisterProcess(msc, particle);
166 
167 
168  // *** Excitation ***
169  ph->RegisterProcess(new G4DNAExcitation("e-_G4DNAExcitation"), particle);
170 
171  // *** Ionisation ***
172  ph->RegisterProcess(new G4DNAIonisation("e-_G4DNAIonisation"), particle);
173 
174  // *** Vibrational excitation ***
175  ph->RegisterProcess(new G4DNAVibExcitation("e-_G4DNAVibExcitation"), particle);
176 
177  // *** Attachment ***
178  ph->RegisterProcess(new G4DNAAttachment("e-_G4DNAAttachment"), particle);
179 
180  } else if ( particleName == "proton" ) {
181 
183  msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
184  ph->RegisterProcess(msc, particle);
185 
186  ph->RegisterProcess(new G4DNAExcitation("proton_G4DNAExcitation"), particle);
187  ph->RegisterProcess(new G4DNAIonisation("proton_G4DNAIonisation"), particle);
188  ph->RegisterProcess(new G4DNAChargeDecrease("proton_G4DNAChargeDecrease"), particle);
189 
190  } else if ( particleName == "hydrogen" ) {
191  ph->RegisterProcess(new G4DNAExcitation("hydrogen_G4DNAExcitation"), particle);
192  ph->RegisterProcess(new G4DNAIonisation("hydrogen_G4DNAIonisation"), particle);
193  ph->RegisterProcess(new G4DNAChargeIncrease("hydrogen_G4DNAChargeIncrease"), particle);
194 
195  } else if ( particleName == "alpha" ) {
196 
198  msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
199  ph->RegisterProcess(msc, particle);
200 
201  ph->RegisterProcess(new G4DNAExcitation("alpha_G4DNAExcitation"), particle);
202  ph->RegisterProcess(new G4DNAIonisation("alpha_G4DNAIonisation"), particle);
203  ph->RegisterProcess(new G4DNAChargeDecrease("alpha_G4DNAChargeDecrease"), particle);
204 
205  } else if ( particleName == "alpha+" ) {
206 
208  msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
209  ph->RegisterProcess(msc, particle);
210 
211  ph->RegisterProcess(new G4DNAExcitation("alpha+_G4DNAExcitation"), particle);
212  ph->RegisterProcess(new G4DNAIonisation("alpha+_G4DNAIonisation"), particle);
213  ph->RegisterProcess(new G4DNAChargeDecrease("alpha+_G4DNAChargeDecrease"), particle);
214  ph->RegisterProcess(new G4DNAChargeIncrease("alpha+_G4DNAChargeIncrease"), particle);
215 
216  } else if ( particleName == "helium" ) {
217  ph->RegisterProcess(new G4DNAExcitation("helium_G4DNAExcitation"), particle);
218  ph->RegisterProcess(new G4DNAIonisation("helium_G4DNAIonisation"), particle);
219  ph->RegisterProcess(new G4DNAChargeIncrease("helium_G4DNAChargeIncrease"), particle);
220 
221  // Extension to HZE proposed by Z. Francis
222 
223  } else if ( particleName == "carbon" ) {
224 
226  msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
227  ph->RegisterProcess(msc, particle);
228 
229  ph->RegisterProcess(new G4DNAIonisation("carbon_G4DNAIonisation"), particle);
230 
231  } else if ( particleName == "nitrogen" ) {
232 
234  msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
235  ph->RegisterProcess(msc, particle);
236 
237  ph->RegisterProcess(new G4DNAIonisation("nitrogen_G4DNAIonisation"), particle);
238 
239  } else if ( particleName == "oxygen" ) {
240 
242  msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
243  ph->RegisterProcess(msc, particle);
244 
245  ph->RegisterProcess(new G4DNAIonisation("oxygen_G4DNAIonisation"), particle);
246 
247  } else if ( particleName == "iron" ) {
248 
250  msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
251  ph->RegisterProcess(msc, particle);
252 
253  ph->RegisterProcess(new G4DNAIonisation("iron_G4DNAIonisation"), particle);
254 
255  }
256 
257  // Warning : the following particles and processes are needed by EM Physics builders
258  // They are taken from the default Livermore Physics list
259  // These particles are currently not handled by Geant4-DNA
260 
261  // e+
262 
263  else if (particleName == "e+") {
264 
265  // Identical to G4EmStandardPhysics_option3
266 
269  G4eIonisation* eIoni = new G4eIonisation();
270  eIoni->SetStepFunction(0.2, 100*um);
271 
272  ph->RegisterProcess(msc, particle);
273  ph->RegisterProcess(eIoni, particle);
274  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
275  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
276 
277  } else if (particleName == "gamma") {
278 
279  G4double LivermoreHighEnergyLimit = GeV;
280 
281  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
282  G4LivermorePhotoElectricModel* theLivermorePhotoElectricModel =
284  theLivermorePhotoElectricModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
285  thePhotoElectricEffect->AddEmModel(0, theLivermorePhotoElectricModel);
286  ph->RegisterProcess(thePhotoElectricEffect, particle);
287 
288  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
289  G4LivermoreComptonModel* theLivermoreComptonModel =
291  theLivermoreComptonModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
292  theComptonScattering->AddEmModel(0, theLivermoreComptonModel);
293  ph->RegisterProcess(theComptonScattering, particle);
294 
295  G4GammaConversion* theGammaConversion = new G4GammaConversion();
296  G4LivermoreGammaConversionModel* theLivermoreGammaConversionModel =
298  theLivermoreGammaConversionModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
299  theGammaConversion->AddEmModel(0, theLivermoreGammaConversionModel);
300  ph->RegisterProcess(theGammaConversion, particle);
301 
302  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
303  G4LivermoreRayleighModel* theRayleighModel = new G4LivermoreRayleighModel();
304  theRayleighModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
305  theRayleigh->AddEmModel(0, theRayleighModel);
306  ph->RegisterProcess(theRayleigh, particle);
307  }
308 
309  // Warning : end of particles and processes are needed by EM Physics builders
310 
311  }
312 
313  // Deexcitation
314  //
317  de->SetFluo(true);
318 }
319 
320 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
static G4GenericIon * GenericIonDefinition()
Definition: G4GenericIon.cc:88
static G4LossTableManager * Instance()
void SetStepFunction(G4double v1, G4double v2)
void SetEmModel(G4VMscModel *, G4int index=1)
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:706
G4GLOB_DLL std::ostream G4cout
#define aParticleIterator
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4Proton * Proton()
Definition: G4Proton.cc:93
static const double GeV
Definition: G4SIunits.hh:196
const G4String & GetPhysicsName() const
static G4DNAGenericIonsManager * Instance(void)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4Positron * Positron()
Definition: G4Positron.cc:94
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=0)
static G4PhysicsListHelper * GetPhysicsListHelper()
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void SetAtomDeexcitation(G4VAtomDeexcitation *)
void SetStepLimitType(G4MscStepLimitType val)
G4ParticleDefinition * GetIon(const G4String &name)
G4_DECLARE_PHYSCONSTR_FACTORY(G4EmDNAPhysics_option1)