Geant4  10.02.p01
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 
74 #include "G4EmParameters.hh"
75 // end of warning
76 
77 #include "G4LossTableManager.hh"
78 #include "G4UAtomicDeexcitation.hh"
79 #include "G4PhysicsListHelper.hh"
80 #include "G4BuilderType.hh"
81 
82 // factory
84 //
86 
87 
88 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
89 
91  : G4VPhysicsConstructor("G4EmDNAPhysics_option1"), verbose(ver)
92 {
95 }
96 
97 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
98 
100  : G4VPhysicsConstructor("G4EmDNAPhysics_option1"), verbose(ver)
101 {
104 }
105 
106 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
107 
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  //genericIonsManager->GetIon("carbon");
134  //genericIonsManager->GetIon("nitrogen");
135  //genericIonsManager->GetIon("oxygen");
136  //genericIonsManager->GetIon("iron");
137 
138 }
139 
140 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
141 
143 {
144  if(verbose > 1) {
145  G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
146  }
148 
149  aParticleIterator->reset();
150  while( (*aParticleIterator)() )
151  {
152  G4ParticleDefinition* particle = aParticleIterator->value();
153  G4String particleName = particle->GetParticleName();
154 
155  if (particleName == "e-") {
156 
157  // *** Elastic scattering (two alternative models available) ***
158 
159  //G4DNAElastic* theDNAElasticProcess = new G4DNAElastic("e-_G4DNAElastic");
160  //theDNAElasticProcess->SetEmModel(new G4DNAChampionElasticModel());
161 
162  // or alternative model
163  //theDNAElasticProcess->SetEmModel(new G4DNAScreenedRutherfordElasticModel());
164 
165  //ph->RegisterProcess(theDNAElasticProcess, particle);
166 
168  msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
169  ph->RegisterProcess(msc, particle);
170 
171 
172  // *** Excitation ***
173  ph->RegisterProcess(new G4DNAExcitation("e-_G4DNAExcitation"), particle);
174 
175  // *** Ionisation ***
176  ph->RegisterProcess(new G4DNAIonisation("e-_G4DNAIonisation"), particle);
177 
178  // *** Vibrational excitation ***
179  ph->RegisterProcess(new G4DNAVibExcitation("e-_G4DNAVibExcitation"), particle);
180 
181  // *** Attachment ***
182  ph->RegisterProcess(new G4DNAAttachment("e-_G4DNAAttachment"), particle);
183 
184  } else if ( particleName == "proton" ) {
185 
187  msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
188  ph->RegisterProcess(msc, particle);
189 
190  ph->RegisterProcess(new G4DNAExcitation("proton_G4DNAExcitation"), particle);
191  ph->RegisterProcess(new G4DNAIonisation("proton_G4DNAIonisation"), particle);
192  ph->RegisterProcess(new G4DNAChargeDecrease("proton_G4DNAChargeDecrease"), particle);
193 
194  } else if ( particleName == "hydrogen" ) {
195  ph->RegisterProcess(new G4DNAExcitation("hydrogen_G4DNAExcitation"), particle);
196  ph->RegisterProcess(new G4DNAIonisation("hydrogen_G4DNAIonisation"), particle);
197  ph->RegisterProcess(new G4DNAChargeIncrease("hydrogen_G4DNAChargeIncrease"), particle);
198 
199  } else if ( particleName == "alpha" ) {
200 
202  msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
203  ph->RegisterProcess(msc, particle);
204 
205  ph->RegisterProcess(new G4DNAExcitation("alpha_G4DNAExcitation"), particle);
206  ph->RegisterProcess(new G4DNAIonisation("alpha_G4DNAIonisation"), particle);
207  ph->RegisterProcess(new G4DNAChargeDecrease("alpha_G4DNAChargeDecrease"), particle);
208 
209  } else if ( particleName == "alpha+" ) {
210 
212  msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
213  ph->RegisterProcess(msc, particle);
214 
215  ph->RegisterProcess(new G4DNAExcitation("alpha+_G4DNAExcitation"), particle);
216  ph->RegisterProcess(new G4DNAIonisation("alpha+_G4DNAIonisation"), particle);
217  ph->RegisterProcess(new G4DNAChargeDecrease("alpha+_G4DNAChargeDecrease"), particle);
218  ph->RegisterProcess(new G4DNAChargeIncrease("alpha+_G4DNAChargeIncrease"), particle);
219 
220  } else if ( particleName == "helium" ) {
221  ph->RegisterProcess(new G4DNAExcitation("helium_G4DNAExcitation"), particle);
222  ph->RegisterProcess(new G4DNAIonisation("helium_G4DNAIonisation"), particle);
223  ph->RegisterProcess(new G4DNAChargeIncrease("helium_G4DNAChargeIncrease"), particle);
224 
225  // Extension to HZE proposed by Z. Francis
226 
227  // S. Incerti
228 
229  } else if ( particleName == "GenericIon" ) {
230 
232  msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
233  ph->RegisterProcess(msc, particle);
234 
235  ph->RegisterProcess(new G4DNAIonisation("GenericIon_G4DNAIonisation"), particle);
236 
237  /*
238  } else if ( particleName == "carbon" ) {
239 
240  G4hMultipleScattering* msc = new G4hMultipleScattering();
241  msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
242  ph->RegisterProcess(msc, particle);
243 
244  ph->RegisterProcess(new G4DNAIonisation("carbon_G4DNAIonisation"), particle);
245 
246  } else if ( particleName == "nitrogen" ) {
247 
248  G4hMultipleScattering* msc = new G4hMultipleScattering();
249  msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
250  ph->RegisterProcess(msc, particle);
251 
252  ph->RegisterProcess(new G4DNAIonisation("nitrogen_G4DNAIonisation"), particle);
253 
254  } else if ( particleName == "oxygen" ) {
255 
256  G4hMultipleScattering* msc = new G4hMultipleScattering();
257  msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
258  ph->RegisterProcess(msc, particle);
259 
260  ph->RegisterProcess(new G4DNAIonisation("oxygen_G4DNAIonisation"), particle);
261 
262  } else if ( particleName == "iron" ) {
263 
264  G4hMultipleScattering* msc = new G4hMultipleScattering();
265  msc->SetEmModel(new G4LowEWentzelVIModel(), 1);
266  ph->RegisterProcess(msc, particle);
267 
268  ph->RegisterProcess(new G4DNAIonisation("iron_G4DNAIonisation"), particle);
269  */
270  }
271 
272  // Warning : the following particles and processes are needed by EM Physics builders
273  // They are taken from the default Livermore Physics list
274  // These particles are currently not handled by Geant4-DNA
275 
276  // e+
277 
278  else if (particleName == "e+") {
279 
280  // Identical to G4EmStandardPhysics_option3
281 
284  G4eIonisation* eIoni = new G4eIonisation();
285  eIoni->SetStepFunction(0.2, 100*um);
286 
287  ph->RegisterProcess(msc, particle);
288  ph->RegisterProcess(eIoni, particle);
289  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
290  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
291 
292  } else if (particleName == "gamma") {
293 
294  G4double LivermoreHighEnergyLimit = GeV;
295 
296  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
297  G4LivermorePhotoElectricModel* theLivermorePhotoElectricModel =
299  theLivermorePhotoElectricModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
300  thePhotoElectricEffect->AddEmModel(0, theLivermorePhotoElectricModel);
301  ph->RegisterProcess(thePhotoElectricEffect, particle);
302 
303  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
304  G4LivermoreComptonModel* theLivermoreComptonModel =
306  theLivermoreComptonModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
307  theComptonScattering->AddEmModel(0, theLivermoreComptonModel);
308  ph->RegisterProcess(theComptonScattering, particle);
309 
310  G4GammaConversion* theGammaConversion = new G4GammaConversion();
311  G4LivermoreGammaConversionModel* theLivermoreGammaConversionModel =
313  theLivermoreGammaConversionModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
314  theGammaConversion->AddEmModel(0, theLivermoreGammaConversionModel);
315  ph->RegisterProcess(theGammaConversion, particle);
316 
317  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
318  G4LivermoreRayleighModel* theRayleighModel = new G4LivermoreRayleighModel();
319  theRayleighModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
320  theRayleigh->AddEmModel(0, theRayleighModel);
321  ph->RegisterProcess(theRayleigh, particle);
322  }
323 
324  // Warning : end of particles and processes are needed by EM Physics builders
325 
326  }
327 
328  // Deexcitation
329  //
332  de->SetFluo(true);
333 }
334 
335 //....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:725
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
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 G4EmParameters * Instance()
static const double um
Definition: G4SIunits.hh:112
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)