Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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: G4EmDNAPhysics.cc 99938 2016-10-12 08:06:52Z gcosmo $
27 // add elastic scattering processes of proton, hydrogen, helium, alpha+, alpha++
28 
29 #include "G4EmDNAPhysics.hh"
30 
31 #include "G4SystemOfUnits.hh"
32 
34 
35 // *** Processes and models for Geant4-DNA
36 
38 #include "G4DNAElastic.hh"
41 #include "G4DNAIonElasticModel.hh"
42 
43 #include "G4DNAExcitation.hh"
44 #include "G4DNAAttachment.hh"
45 #include "G4DNAVibExcitation.hh"
46 #include "G4DNAIonisation.hh"
47 #include "G4DNAChargeDecrease.hh"
48 #include "G4DNAChargeIncrease.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 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
88 
90  : G4VPhysicsConstructor("G4EmDNAPhysics"), verbose(ver)
91 {
93  param->SetDefaults();
94  param->SetFluo(true);
95  param->SetAuger(true);
96  param->SetAugerCascade(true);
97  param->SetDeexcitationIgnoreCut(true);
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  auto myParticleIterator=GetParticleIterator();
146  myParticleIterator->reset();
147  while( (*myParticleIterator)() )
148  {
149  G4ParticleDefinition* particle = myParticleIterator->value();
150  G4String particleName = particle->GetParticleName();
151 
152  if (particleName == "e-") {
153 
154  G4DNAElectronSolvation* solvation =
155  new G4DNAElectronSolvation("e-_G4DNAElectronSolvation");
158  therm->SetHighEnergyLimit(7.4*eV); // limit of the Champion's model
159  solvation->SetEmModel(therm);
160  ph->RegisterProcess(solvation, particle);
161 
162  // *** Elastic scattering (two alternative models available) ***
163 
164  G4DNAElastic* theDNAElasticProcess = new G4DNAElastic("e-_G4DNAElastic");
165  theDNAElasticProcess->SetEmModel(new G4DNAChampionElasticModel());
166 
167  // or alternative model
168  //theDNAElasticProcess->SetEmModel(new G4DNAScreenedRutherfordElasticModel());
169 
170  ph->RegisterProcess(theDNAElasticProcess, particle);
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  //ph->RegisterProcess(new G4DNAElastic("proton_G4DNAElastic"), particle);
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 G4DNAElastic("hydrogen_G4DNAElastic"), particle);
192  ph->RegisterProcess(new G4DNAExcitation("hydrogen_G4DNAExcitation"), particle);
193  ph->RegisterProcess(new G4DNAIonisation("hydrogen_G4DNAIonisation"), particle);
194  ph->RegisterProcess(new G4DNAChargeIncrease("hydrogen_G4DNAChargeIncrease"), particle);
195 
196  } else if ( particleName == "alpha" ) {
197  //ph->RegisterProcess(new G4DNAElastic("alpha_G4DNAElastic"), particle);
198  ph->RegisterProcess(new G4DNAExcitation("alpha_G4DNAExcitation"), particle);
199  ph->RegisterProcess(new G4DNAIonisation("alpha_G4DNAIonisation"), particle);
200  ph->RegisterProcess(new G4DNAChargeDecrease("alpha_G4DNAChargeDecrease"), particle);
201 
202  } else if ( particleName == "alpha+" ) {
203  //ph->RegisterProcess(new G4DNAElastic("alpha+_G4DNAElastic"), particle);
204  ph->RegisterProcess(new G4DNAExcitation("alpha+_G4DNAExcitation"), particle);
205  ph->RegisterProcess(new G4DNAIonisation("alpha+_G4DNAIonisation"), particle);
206  ph->RegisterProcess(new G4DNAChargeDecrease("alpha+_G4DNAChargeDecrease"), particle);
207  ph->RegisterProcess(new G4DNAChargeIncrease("alpha+_G4DNAChargeIncrease"), particle);
208 
209  } else if ( particleName == "helium" ) {
210  //ph->RegisterProcess(new G4DNAElastic("helium_G4DNAElastic"), particle);
211  ph->RegisterProcess(new G4DNAExcitation("helium_G4DNAExcitation"), particle);
212  ph->RegisterProcess(new G4DNAIonisation("helium_G4DNAIonisation"), particle);
213  ph->RegisterProcess(new G4DNAChargeIncrease("helium_G4DNAChargeIncrease"), particle);
214 
215  // Extension to HZE proposed by Z. Francis
216 
217  //SEB
218  } else if ( particleName == "GenericIon" ) {
219  ph->RegisterProcess(new G4DNAIonisation("GenericIon_G4DNAIonisation"), particle);
220 
221  /*
222  } else if ( particleName == "carbon" ) {
223  ph->RegisterProcess(new G4DNAIonisation("carbon_G4DNAIonisation"), particle);
224 
225  } else if ( particleName == "nitrogen" ) {
226  ph->RegisterProcess(new G4DNAIonisation("nitrogen_G4DNAIonisation"), particle);
227 
228  } else if ( particleName == "oxygen" ) {
229  ph->RegisterProcess(new G4DNAIonisation("oxygen_G4DNAIonisation"), particle);
230 
231  } else if ( particleName == "iron" ) {
232  ph->RegisterProcess(new G4DNAIonisation("iron_G4DNAIonisation"), particle);
233  */
234 
235  }
236 
237  // Warning : the following particles and processes are needed by EM Physics builders
238  // They are taken from the default Livermore Physics list
239  // These particles are currently not handled by Geant4-DNA
240 
241  // e+
242 
243  else if (particleName == "e+") {
244 
245  // Identical to G4EmStandardPhysics_option3
246 
249  G4eIonisation* eIoni = new G4eIonisation();
250  eIoni->SetStepFunction(0.2, 100*um);
251 
252  ph->RegisterProcess(msc, particle);
253  ph->RegisterProcess(eIoni, particle);
254  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
255  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
256 
257  } else if (particleName == "gamma") {
258 
259  G4double LivermoreHighEnergyLimit = GeV;
260 
261  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
262  G4LivermorePhotoElectricModel* theLivermorePhotoElectricModel =
264  theLivermorePhotoElectricModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
265  thePhotoElectricEffect->AddEmModel(0, theLivermorePhotoElectricModel);
266  ph->RegisterProcess(thePhotoElectricEffect, particle);
267 
268  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
269  G4LivermoreComptonModel* theLivermoreComptonModel =
271  theLivermoreComptonModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
272  theComptonScattering->AddEmModel(0, theLivermoreComptonModel);
273  ph->RegisterProcess(theComptonScattering, particle);
274 
275  G4GammaConversion* theGammaConversion = new G4GammaConversion();
276  G4LivermoreGammaConversionModel* theLivermoreGammaConversionModel =
278  theLivermoreGammaConversionModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
279  theGammaConversion->AddEmModel(0, theLivermoreGammaConversionModel);
280  ph->RegisterProcess(theGammaConversion, particle);
281 
282  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
283  G4LivermoreRayleighModel* theRayleighModel = new G4LivermoreRayleighModel();
284  theRayleighModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
285  theRayleigh->AddEmModel(0, theRayleighModel);
286  ph->RegisterProcess(theRayleigh, particle);
287  }
288 
289  // Warning : end of particles and processes are needed by EM Physics builders
290 
291  }
292 
293  // Deexcitation
294  //
297 }
298 
299 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
virtual void ConstructParticle()
static G4GenericIon * GenericIonDefinition()
Definition: G4GenericIon.cc:88
static G4LossTableManager * Instance()
void SetDeexcitationIgnoreCut(G4bool val)
virtual ~G4EmDNAPhysics()
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:732
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
void SetEmModel(G4VEmModel *, G4int index=1)
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
virtual void ConstructProcess()
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
G4EmDNAPhysics(G4int ver=1, const G4String &name="")
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)