Geant4  10.02.p01
G4EmDNAPhysics_option7.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 // S. Incerti (incerti@cenbg.in2p3.fr)
27 //
28 
30 
31 #include "G4SystemOfUnits.hh"
32 
34 
35 // *** Processes and models for Geant4-DNA
36 
37 #include "G4DNAElastic.hh"
41 
42 #include "G4DNAIonisation.hh"
44 
45 #include "G4DNAExcitation.hh"
47 
48 #include "G4DNAAttachment.hh"
49 #include "G4DNAVibExcitation.hh"
50 
51 #include "G4DNAChargeDecrease.hh"
52 #include "G4DNAChargeIncrease.hh"
53 
54 // particles
55 
56 #include "G4Electron.hh"
57 #include "G4Proton.hh"
58 #include "G4GenericIon.hh"
59 
60 // Warning : the following is needed in order to use EM Physics builders
61 // e+
62 #include "G4Positron.hh"
63 #include "G4eMultipleScattering.hh"
64 #include "G4eIonisation.hh"
65 #include "G4eBremsstrahlung.hh"
66 #include "G4eplusAnnihilation.hh"
67 // gamma
68 #include "G4Gamma.hh"
69 #include "G4PhotoElectricEffect.hh"
71 #include "G4ComptonScattering.hh"
73 #include "G4GammaConversion.hh"
75 #include "G4RayleighScattering.hh"
77 // end of warning
78 
79 #include "G4LossTableManager.hh"
80 #include "G4UAtomicDeexcitation.hh"
81 #include "G4PhysicsListHelper.hh"
82 #include "G4BuilderType.hh"
83 
84 // factory
86 //
88 
89 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
90 
92  G4VPhysicsConstructor("G4EmDNAPhysics_option7"), verbose(ver)
93 {
95 }
96 
97 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
98 
100  G4VPhysicsConstructor("G4EmDNAPhysics_option7"), verbose(ver)
101 {
103 }
104 
105 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
106 
108 {
109 }
110 
111 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
112 
114 {
115 // bosons
116  G4Gamma::Gamma();
117 
118 // leptons
121 
122 // baryons
124 
126 
127  G4DNAGenericIonsManager * genericIonsManager;
128  genericIonsManager = G4DNAGenericIonsManager::Instance();
129  genericIonsManager->GetIon("alpha++");
130  genericIonsManager->GetIon("alpha+");
131  genericIonsManager->GetIon("helium");
132  genericIonsManager->GetIon("hydrogen");
133 
134 }
135 
136 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
137 
139 {
140  if(verbose > 1)
141  {
142  G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
143  }
145 
146  aParticleIterator->reset();
147  while( (*aParticleIterator)() )
148  {
149  G4ParticleDefinition* particle = aParticleIterator->value();
150  G4String particleName = particle->GetParticleName();
151 
152  if (particleName == "e-")
153  {
154 
155  // *** Elastic scattering (two alternative models available) ***
156  G4DNAElastic* theDNAElasticProcess = new G4DNAElastic("e-_G4DNAElastic");
157  theDNAElasticProcess->SetEmModel(new G4DNAUeharaScreenedRutherfordElasticModel());
158 
159  ph->RegisterProcess(theDNAElasticProcess, particle);
160 
161  {
162  // *** Excitation ***
163  G4DNAExcitation* theDNAExcitationProcess =
164  new G4DNAExcitation("e-_G4DNAExcitation");
165  {
167  bornExc->SetLowEnergyLimit(10*keV);
168  bornExc->SetHighEnergyLimit(1.*MeV);
169  bornExc->SetActivationLowEnergyLimit(10*keV);
170  bornExc->SetActivationHighEnergyLimit(1.*MeV);
171  theDNAExcitationProcess->SetEmModel(bornExc,1);
172  theDNAExcitationProcess->AddEmModel(1, bornExc);
173  }
174  {
176  emExc->SetActivationLowEnergyLimit(0);
177  emExc->SetActivationHighEnergyLimit(10.*keV);
178  theDNAExcitationProcess->SetEmModel(emExc,2);
179  theDNAExcitationProcess->AddEmModel(2, emExc);
180  }
181  ph->RegisterProcess(theDNAExcitationProcess, particle);
182  }
183  {
184  // *** Ionisation ***
185  G4DNAIonisation* theDNAIonisationProcess =
186  new G4DNAIonisation("e-_G4DNAIonisation");
187  {
189  bornIon->SetLowEnergyLimit(10*keV);
190  bornIon->SetHighEnergyLimit(1.*MeV);
191  bornIon->SetActivationLowEnergyLimit(10*keV);
192  bornIon->SetActivationHighEnergyLimit(1.*MeV);
193 
194  theDNAIonisationProcess->SetEmModel(bornIon,1);
195  theDNAIonisationProcess->AddEmModel(1,bornIon);
196  }
197 
198  {
200  emIonModel->SetLowEnergyLimit(10*eV);
201  emIonModel->SetHighEnergyLimit(10.*keV);
202  emIonModel->SetActivationLowEnergyLimit(0);
203  emIonModel->SetActivationHighEnergyLimit(10.*keV);
204 
205  theDNAIonisationProcess->SetEmModel(emIonModel,2);
206  emIonModel->SelectFasterComputation(true);
207  theDNAIonisationProcess->AddEmModel(2,emIonModel);
208  }
209  ph->RegisterProcess(theDNAIonisationProcess, particle);
210  }
211  // *** Vibrational excitation ***
212  ph->RegisterProcess(new G4DNAVibExcitation("e-_G4DNAVibExcitation"), particle);
213 
214  // *** Attachment ***
215  ph->RegisterProcess(new G4DNAAttachment("e-_G4DNAAttachment"), particle);
216 
217  } else if ( particleName == "proton" ) {
218  ph->RegisterProcess(new G4DNAElastic("proton_G4DNAElastic"), particle);
219  ph->RegisterProcess(new G4DNAExcitation("proton_G4DNAExcitation"), particle);
220  ph->RegisterProcess(new G4DNAIonisation("proton_G4DNAIonisation"), particle);
221  ph->RegisterProcess(new G4DNAChargeDecrease("proton_G4DNAChargeDecrease"), particle);
222 
223  } else if ( particleName == "hydrogen" ) {
224  ph->RegisterProcess(new G4DNAElastic("hydrogen_G4DNAElastic"), particle);
225  ph->RegisterProcess(new G4DNAExcitation("hydrogen_G4DNAExcitation"), particle);
226  ph->RegisterProcess(new G4DNAIonisation("hydrogen_G4DNAIonisation"), particle);
227  ph->RegisterProcess(new G4DNAChargeIncrease("hydrogen_G4DNAChargeIncrease"), particle);
228 
229  } else if ( particleName == "alpha" ) {
230  ph->RegisterProcess(new G4DNAElastic("alpha_G4DNAElastic"), particle);
231  ph->RegisterProcess(new G4DNAExcitation("alpha_G4DNAExcitation"), particle);
232  ph->RegisterProcess(new G4DNAIonisation("alpha_G4DNAIonisation"), particle);
233  ph->RegisterProcess(new G4DNAChargeDecrease("alpha_G4DNAChargeDecrease"), particle);
234 
235  } else if ( particleName == "alpha+" ) {
236  ph->RegisterProcess(new G4DNAElastic("alpha+_G4DNAElastic"), particle);
237  ph->RegisterProcess(new G4DNAExcitation("alpha+_G4DNAExcitation"), particle);
238  ph->RegisterProcess(new G4DNAIonisation("alpha+_G4DNAIonisation"), particle);
239  ph->RegisterProcess(new G4DNAChargeDecrease("alpha+_G4DNAChargeDecrease"), particle);
240  ph->RegisterProcess(new G4DNAChargeIncrease("alpha+_G4DNAChargeIncrease"), particle);
241 
242  } else if ( particleName == "helium" ) {
243  ph->RegisterProcess(new G4DNAElastic("helium_G4DNAElastic"), particle);
244  ph->RegisterProcess(new G4DNAExcitation("helium_G4DNAExcitation"), particle);
245  ph->RegisterProcess(new G4DNAIonisation("helium_G4DNAIonisation"), particle);
246  ph->RegisterProcess(new G4DNAChargeIncrease("helium_G4DNAChargeIncrease"), particle);
247  }
248  // Extension to HZE proposed by Z. Francis
249  else if ( particleName == "GenericIon" ) {
250  ph->RegisterProcess(new G4DNAIonisation("GenericIon_G4DNAIonisation"), particle);
251  }
252 
253  // Warning : the following particles and processes are needed by EM Physics builders
254  // They are taken from the default Livermore Physics list
255  // These particles are currently not handled by Geant4-DNA
256 
257  // e+
258 
259  else if (particleName == "e+") {
260 
261  // Identical to G4EmStandardPhysics_option3
262 
265  G4eIonisation* eIoni = new G4eIonisation();
266  eIoni->SetStepFunction(0.2, 100*um);
267 
268  ph->RegisterProcess(msc, particle);
269  ph->RegisterProcess(eIoni, particle);
270  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
271  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
272 
273  } else if (particleName == "gamma") {
274 
275  G4double LivermoreHighEnergyLimit = GeV;
276 
277  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
278  G4LivermorePhotoElectricModel* theLivermorePhotoElectricModel =
280  theLivermorePhotoElectricModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
281  thePhotoElectricEffect->AddEmModel(0, theLivermorePhotoElectricModel);
282  ph->RegisterProcess(thePhotoElectricEffect, particle);
283 
284  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
285  G4LivermoreComptonModel* theLivermoreComptonModel =
287  theLivermoreComptonModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
288  theComptonScattering->AddEmModel(0, theLivermoreComptonModel);
289  ph->RegisterProcess(theComptonScattering, particle);
290 
291  G4GammaConversion* theGammaConversion = new G4GammaConversion();
292  G4LivermoreGammaConversionModel* theLivermoreGammaConversionModel =
294  theLivermoreGammaConversionModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
295  theGammaConversion->AddEmModel(0, theLivermoreGammaConversionModel);
296  ph->RegisterProcess(theGammaConversion, particle);
297 
298  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
299  G4LivermoreRayleighModel* theRayleighModel = new G4LivermoreRayleighModel();
300  theRayleighModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
301  theRayleigh->AddEmModel(0, theRayleighModel);
302  ph->RegisterProcess(theRayleigh, particle);
303  }
304 
305  // Warning : end of particles and processes are needed by EM Physics builders
306 
307  }
308 
309  // Deexcitation
310  //
313  de->SetFluo(true);
314 }
315 
316 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
void SetActivationHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:739
static G4GenericIon * GenericIonDefinition()
Definition: G4GenericIon.cc:88
static const double MeV
Definition: G4SIunits.hh:211
#define G4DNABornIonisationModel
static G4LossTableManager * Instance()
void SetStepFunction(G4double v1, G4double v2)
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4DNABornExcitationModel1 G4DNABornExcitationModel
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
void SetEmModel(G4VEmModel *, G4int index=1)
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:214
G4_DECLARE_PHYSCONSTR_FACTORY(G4EmDNAPhysics_option7)
const G4String & GetPhysicsName() const
static G4DNAGenericIonsManager * Instance(void)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
static const double eV
Definition: G4SIunits.hh:212
static G4Positron * Positron()
Definition: G4Positron.cc:94
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=0)
static G4PhysicsListHelper * GetPhysicsListHelper()
static const double um
Definition: G4SIunits.hh:112
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
static const double keV
Definition: G4SIunits.hh:213
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
void SetAtomDeexcitation(G4VAtomDeexcitation *)
void SetStepLimitType(G4MscStepLimitType val)
G4ParticleDefinition * GetIon(const G4String &name)