Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4EmDNAPhysics_option3.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 // added elastic scattering processes of proton, hydrogen, helium, alpha+, alpha
27 
29 
30 #include "G4SystemOfUnits.hh"
31 
33 
34 // *** Processes and models for Geant4-DNA
35 
37 #include "G4DNAElastic.hh"
40 #include "G4DNAIonElasticModel.hh"
41 
42 #include "G4DNAExcitation.hh"
43 #include "G4DNAAttachment.hh"
44 #include "G4DNAVibExcitation.hh"
45 #include "G4DNAIonisation.hh"
46 #include "G4DNAChargeDecrease.hh"
47 #include "G4DNAChargeIncrease.hh"
48 
49 // particles
50 
51 #include "G4Electron.hh"
52 #include "G4Proton.hh"
53 #include "G4GenericIon.hh"
54 
55 // Warning : the following is needed in order to use EM Physics builders
56 // e+
57 #include "G4Positron.hh"
58 #include "G4eMultipleScattering.hh"
59 #include "G4eIonisation.hh"
60 #include "G4eBremsstrahlung.hh"
61 #include "G4eplusAnnihilation.hh"
62 // gamma
63 #include "G4Gamma.hh"
64 #include "G4PhotoElectricEffect.hh"
66 #include "G4ComptonScattering.hh"
68 #include "G4GammaConversion.hh"
70 #include "G4RayleighScattering.hh"
72 
73 #include "G4EmParameters.hh"
74 // end of warning
75 
76 #include "G4LossTableManager.hh"
77 #include "G4UAtomicDeexcitation.hh"
78 #include "G4PhysicsListHelper.hh"
79 #include "G4BuilderType.hh"
80 
81 // factory
83 //
85 
86 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
87 
89  : G4VPhysicsConstructor("G4EmDNAPhysics_option3"), verbose(ver)
90 {
92  param->SetDefaults();
93  param->SetFluo(true);
94  param->SetAuger(true);
95  param->SetAugerCascade(true);
96  param->SetDeexcitationIgnoreCut(true);
97 
99 }
100 
101 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
102 
104 {}
105 
106 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
107 
109 {
110 // bosons
111  G4Gamma::Gamma();
112 
113 // leptons
116 
117 // baryons
119 
121 
122  G4DNAGenericIonsManager * genericIonsManager;
123  genericIonsManager=G4DNAGenericIonsManager::Instance();
124  genericIonsManager->GetIon("alpha++");
125  genericIonsManager->GetIon("alpha+");
126  genericIonsManager->GetIon("helium");
127  genericIonsManager->GetIon("hydrogen");
128 
129 }
130 
131 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
132 
134 {
135  if(verbose > 1) {
136  G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
137  }
139 
140  auto myParticleIterator=GetParticleIterator();
141  myParticleIterator->reset();
142  while( (*myParticleIterator)() )
143  {
144  G4ParticleDefinition* particle = myParticleIterator->value();
145  G4String particleName = particle->GetParticleName();
146 
147  if (particleName == "e-") {
148 
149  G4DNAElectronSolvation* solvation =
150  new G4DNAElectronSolvation("e-_G4DNAElectronSolvation");
153  therm->SetHighEnergyLimit(7.4*eV); // limit of the Champion's model
154  solvation->SetEmModel(therm);
155  ph->RegisterProcess(solvation, particle);
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 
167  // *** Excitation ***
168  ph->RegisterProcess(new G4DNAExcitation("e-_G4DNAExcitation"), particle);
169 
170  // *** Ionisation ***
171  ph->RegisterProcess(new G4DNAIonisation("e-_G4DNAIonisation"), particle);
172 
173  // *** Vibrational excitation ***
174  ph->RegisterProcess(new G4DNAVibExcitation("e-_G4DNAVibExcitation"), particle);
175 
176  // *** Attachment ***
177  ph->RegisterProcess(new G4DNAAttachment("e-_G4DNAAttachment"), particle);
178 
179  } else if ( particleName == "proton" ) {
180  ph->RegisterProcess(new G4DNAElastic("proton_G4DNAElastic"), particle);
181  ph->RegisterProcess(new G4DNAExcitation("proton_G4DNAExcitation"), particle);
182  ph->RegisterProcess(new G4DNAIonisation("proton_G4DNAIonisation"), particle);
183  ph->RegisterProcess(new G4DNAChargeDecrease("proton_G4DNAChargeDecrease"), particle);
184 
185  } else if ( particleName == "hydrogen" ) {
186  ph->RegisterProcess(new G4DNAElastic("hydrogen_G4DNAElastic"), particle);
187  ph->RegisterProcess(new G4DNAExcitation("hydrogen_G4DNAExcitation"), particle);
188  ph->RegisterProcess(new G4DNAIonisation("hydrogen_G4DNAIonisation"), particle);
189  ph->RegisterProcess(new G4DNAChargeIncrease("hydrogen_G4DNAChargeIncrease"), particle);
190 
191  } else if ( particleName == "alpha" ) {
192  ph->RegisterProcess(new G4DNAElastic("alpha_G4DNAElastic"), particle);
193  ph->RegisterProcess(new G4DNAExcitation("alpha_G4DNAExcitation"), particle);
194  ph->RegisterProcess(new G4DNAIonisation("alpha_G4DNAIonisation"), particle);
195  ph->RegisterProcess(new G4DNAChargeDecrease("alpha_G4DNAChargeDecrease"), particle);
196 
197  } else if ( particleName == "alpha+" ) {
198  ph->RegisterProcess(new G4DNAElastic("alpha+_G4DNAElastic"), particle);
199  ph->RegisterProcess(new G4DNAExcitation("alpha+_G4DNAExcitation"), particle);
200  ph->RegisterProcess(new G4DNAIonisation("alpha+_G4DNAIonisation"), particle);
201  ph->RegisterProcess(new G4DNAChargeDecrease("alpha+_G4DNAChargeDecrease"), particle);
202  ph->RegisterProcess(new G4DNAChargeIncrease("alpha+_G4DNAChargeIncrease"), particle);
203 
204  } else if ( particleName == "helium" ) {
205  ph->RegisterProcess(new G4DNAElastic("helium_G4DNAElastic"), particle);
206  ph->RegisterProcess(new G4DNAExcitation("helium_G4DNAExcitation"), particle);
207  ph->RegisterProcess(new G4DNAIonisation("helium_G4DNAIonisation"), particle);
208  ph->RegisterProcess(new G4DNAChargeIncrease("helium_G4DNAChargeIncrease"), particle);
209 
210  } else if ( particleName == "GenericIon" ) {
211  ph->RegisterProcess(new G4DNAIonisation("GenericIon_G4DNAIonisation"), particle);
212 
213  }
214 
215  // Warning : the following particles and processes are needed by EM Physics builders
216  // They are taken from the default Livermore Physics list
217  // These particles are currently not handled by Geant4-DNA
218 
219  // e+
220 
221  else if (particleName == "e+") {
222 
223  // Identical to G4EmStandardPhysics_option3
224 
227  G4eIonisation* eIoni = new G4eIonisation();
228  eIoni->SetStepFunction(0.2, 100*um);
229 
230  ph->RegisterProcess(msc, particle);
231  ph->RegisterProcess(eIoni, particle);
232  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
233  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
234 
235  } else if (particleName == "gamma") {
236 
237  G4double LivermoreHighEnergyLimit = GeV;
238 
239  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
240  G4LivermorePhotoElectricModel* theLivermorePhotoElectricModel =
242  theLivermorePhotoElectricModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
243  thePhotoElectricEffect->AddEmModel(0, theLivermorePhotoElectricModel);
244  ph->RegisterProcess(thePhotoElectricEffect, particle);
245 
246  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
247  G4LivermoreComptonModel* theLivermoreComptonModel =
249  theLivermoreComptonModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
250  theComptonScattering->AddEmModel(0, theLivermoreComptonModel);
251  ph->RegisterProcess(theComptonScattering, particle);
252 
253  G4GammaConversion* theGammaConversion = new G4GammaConversion();
254  G4LivermoreGammaConversionModel* theLivermoreGammaConversionModel =
256  theLivermoreGammaConversionModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
257  theGammaConversion->AddEmModel(0, theLivermoreGammaConversionModel);
258  ph->RegisterProcess(theGammaConversion, particle);
259 
260  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
261  G4LivermoreRayleighModel* theRayleighModel = new G4LivermoreRayleighModel();
262  theRayleighModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
263  theRayleigh->AddEmModel(0, theRayleighModel);
264  ph->RegisterProcess(theRayleigh, particle);
265  }
266 
267  // Warning : end of particles and processes are needed by EM Physics builders
268 
269  }
270 
271  // Deexcitation
272  //
275 }
276 
277 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
static G4GenericIon * GenericIonDefinition()
Definition: G4GenericIon.cc:88
static G4LossTableManager * Instance()
void SetDeexcitationIgnoreCut(G4bool val)
void SetAuger(G4bool val)
int G4int
Definition: G4Types.hh:78
void SetStepFunction(G4double v1, G4double v2, G4bool lock=true)
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:724
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
void SetEmModel(G4VEmModel *, G4int index=1)
G4EmDNAPhysics_option3(G4int ver=1, const G4String &name="")
G4GLOB_DLL std::ostream G4cout
static constexpr double um
Definition: G4SIunits.hh:113
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
void SetAugerCascade(G4bool val)
static G4Proton * Proton()
Definition: G4Proton.cc:93
static constexpr double eV
Definition: G4SIunits.hh:215
const G4String & GetPhysicsName() const
static G4DNAGenericIonsManager * Instance(void)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
G4TDNAOneStepThermalizationModel< DNA::Penetration::Meesungnoen2002 > G4DNAOneStepThermalizationModel
static G4Positron * Positron()
Definition: G4Positron.cc:94
static G4PhysicsListHelper * GetPhysicsListHelper()
static G4EmParameters * Instance()
static constexpr double GeV
Definition: G4SIunits.hh:217
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void SetAtomDeexcitation(G4VAtomDeexcitation *)
void SetFluo(G4bool val)
void SetStepLimitType(G4MscStepLimitType val)
G4ParticleDefinition * GetIon(const G4String &name)
#define G4_DECLARE_PHYSCONSTR_FACTORY(physics_constructor)