Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4EmDNAPhysics_option2.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_option2.cc 82360 2014-06-17 15:16:30Z gcosmo $
27 
28 // SI: This constructor uses speedup options of DNA models
29 
31 
32 #include "G4SystemOfUnits.hh"
33 
35 
36 // *** Processes and models for Geant4-DNA
37 
39 #include "G4DNAElastic.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_option2"), verbose(ver)
91 {
93 
94  param->SetDefaults();
95  param->SetFluo(true);
96  param->SetAuger(true);
97  param->SetAugerCascade(true);
98  param->SetDeexcitationIgnoreCut(true);
99 
101 }
102 
103 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
104 
106 {}
107 
108 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
109 
111 {
112 // bosons
113  G4Gamma::Gamma();
114 
115 // leptons
118 
119 // baryons
121 
123 
124  G4DNAGenericIonsManager * genericIonsManager;
125  genericIonsManager=G4DNAGenericIonsManager::Instance();
126  genericIonsManager->GetIon("alpha++");
127  genericIonsManager->GetIon("alpha+");
128  genericIonsManager->GetIon("helium");
129  genericIonsManager->GetIon("hydrogen");
130 
131 }
132 
133 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
134 
136 {
137  if(verbose > 1) {
138  G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
139  }
141 
142  auto myParticleIterator=GetParticleIterator();
143  myParticleIterator->reset();
144  while( (*myParticleIterator)() )
145  {
146  G4ParticleDefinition* particle = myParticleIterator->value();
147  G4String particleName = particle->GetParticleName();
148 
149  if (particleName == "e-") {
150 
151  G4DNAElectronSolvation* solvation =
152  new G4DNAElectronSolvation("e-_G4DNAElectronSolvation");
155  therm->SetHighEnergyLimit(7.4*eV); // limit of the Champion's model
156  solvation->SetEmModel(therm);
157  ph->RegisterProcess(solvation, particle);
158 
159  // *** Elastic scattering (two alternative models available) ***
160 
161  G4DNAElastic* theDNAElasticProcess = new G4DNAElastic("e-_G4DNAElastic");
162  theDNAElasticProcess->SetEmModel(new G4DNAChampionElasticModel());
163 
164  // or alternative model
165  //theDNAElasticProcess->SetEmModel(new G4DNAScreenedRutherfordElasticModel());
166 
167  ph->RegisterProcess(theDNAElasticProcess, particle);
168 
169  // *** Excitation ***
170  ph->RegisterProcess(new G4DNAExcitation("e-_G4DNAExcitation"), particle);
171 
172  // *** Ionisation ***
173  //ph->RegisterProcess(new G4DNAIonisation("e-_G4DNAIonisation"), particle);
174  G4DNAIonisation* theDNAIonisationProcess = new G4DNAIonisation("e-_G4DNAIonisation");
175  theDNAIonisationProcess->SetEmModel(new G4DNABornIonisationModel());
176  ((G4DNABornIonisationModel*)(theDNAIonisationProcess->EmModel()))->SelectFasterComputation(true);
177  ph->RegisterProcess(theDNAIonisationProcess, particle);
178 
179  // *** Vibrational excitation ***
180  ph->RegisterProcess(new G4DNAVibExcitation("e-_G4DNAVibExcitation"), particle);
181 
182  // *** Attachment ***
183  ph->RegisterProcess(new G4DNAAttachment("e-_G4DNAAttachment"), particle);
184 
185  } else if ( particleName == "proton" ) {
186 
187  ph->RegisterProcess(new G4DNAExcitation("proton_G4DNAExcitation"), particle);
188 
189  G4DNAIonisation* theDNAIonisationProcess = new G4DNAIonisation("proton_G4DNAIonisation");
190 
191  G4VEmModel* mod1;
193  mod1->SetLowEnergyLimit(0*eV);
194  mod1->SetHighEnergyLimit(500*keV);
195 
196  G4VEmModel* mod2;
197  mod2= new G4DNABornIonisationModel();
198  mod2->SetLowEnergyLimit(500*keV);
199  mod2->SetHighEnergyLimit(100*MeV);
200 
201  theDNAIonisationProcess->SetEmModel(mod1,1);
202  theDNAIonisationProcess->SetEmModel(mod2,2);
203  ((G4DNABornIonisationModel*)(theDNAIonisationProcess->EmModel(2)))->SelectFasterComputation(true);
204 
205  ph->RegisterProcess(theDNAIonisationProcess, particle);
206 
207  ph->RegisterProcess(new G4DNAChargeDecrease("proton_G4DNAChargeDecrease"), particle);
208 
209  } else if ( particleName == "hydrogen" ) {
210 
211  ph->RegisterProcess(new G4DNAExcitation("hydrogen_G4DNAExcitation"), particle);
212 
213  //ph->RegisterProcess(new G4DNAIonisation("hydrogen_G4DNAIonisation"), particle);
214  G4DNAIonisation* theDNAIonisationProcess = new G4DNAIonisation("hydrogen_G4DNAIonisation");
215  theDNAIonisationProcess->SetEmModel(new G4DNARuddIonisationExtendedModel());
216  ph->RegisterProcess(theDNAIonisationProcess, particle);
217 
218  ph->RegisterProcess(new G4DNAChargeIncrease("hydrogen_G4DNAChargeIncrease"), particle);
219 
220  } else if ( particleName == "alpha" ) {
221 
222  ph->RegisterProcess(new G4DNAExcitation("alpha_G4DNAExcitation"), particle);
223 
224  //ph->RegisterProcess(new G4DNAIonisation("alpha_G4DNAIonisation"), particle);
225  G4DNAIonisation* theDNAIonisationProcess = new G4DNAIonisation("alpha_G4DNAIonisation");
226  theDNAIonisationProcess->SetEmModel(new G4DNARuddIonisationExtendedModel());
227  ph->RegisterProcess(theDNAIonisationProcess, particle);
228 
229  ph->RegisterProcess(new G4DNAChargeDecrease("alpha_G4DNAChargeDecrease"), particle);
230 
231  } else if ( particleName == "alpha+" ) {
232 
233  ph->RegisterProcess(new G4DNAExcitation("alpha+_G4DNAExcitation"), particle);
234 
235  //ph->RegisterProcess(new G4DNAIonisation("alpha+_G4DNAIonisation"), particle);
236  G4DNAIonisation* theDNAIonisationProcess = new G4DNAIonisation("alpha+_G4DNAIonisation");
237  theDNAIonisationProcess->SetEmModel(new G4DNARuddIonisationExtendedModel());
238  ph->RegisterProcess(theDNAIonisationProcess, particle);
239 
240  ph->RegisterProcess(new G4DNAChargeDecrease("alpha+_G4DNAChargeDecrease"), particle);
241  ph->RegisterProcess(new G4DNAChargeIncrease("alpha+_G4DNAChargeIncrease"), particle);
242 
243  } else if ( particleName == "helium" ) {
244 
245  ph->RegisterProcess(new G4DNAExcitation("helium_G4DNAExcitation"), particle);
246 
247  //ph->RegisterProcess(new G4DNAIonisation("helium_G4DNAIonisation"), particle);
248  G4DNAIonisation* theDNAIonisationProcess = new G4DNAIonisation("helium_G4DNAIonisation");
249  theDNAIonisationProcess->SetEmModel(new G4DNARuddIonisationExtendedModel());
250  ph->RegisterProcess(theDNAIonisationProcess, particle);
251 
252  ph->RegisterProcess(new G4DNAChargeIncrease("helium_G4DNAChargeIncrease"), particle);
253 
254  } else if ( particleName == "GenericIon" ) {
255  ph->RegisterProcess(new G4DNAIonisation("GenericIon_G4DNAIonisation"), particle);
256 
257  }
258 
259  // Warning : the following particles and processes are needed by EM Physics builders
260  // They are taken from the default Livermore Physics list
261  // These particles are currently not handled by Geant4-DNA
262 
263  // e+
264 
265  else if (particleName == "e+") {
266 
267  // Identical to G4EmStandardPhysics_option3
268 
271  G4eIonisation* eIoni = new G4eIonisation();
272  eIoni->SetStepFunction(0.2, 100*um);
273 
274  ph->RegisterProcess(msc, particle);
275  ph->RegisterProcess(eIoni, particle);
276  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
277  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
278 
279  } else if (particleName == "gamma") {
280 
281  G4double LivermoreHighEnergyLimit = GeV;
282 
283  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
284  G4LivermorePhotoElectricModel* theLivermorePhotoElectricModel =
286  theLivermorePhotoElectricModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
287  thePhotoElectricEffect->AddEmModel(0, theLivermorePhotoElectricModel);
288  ph->RegisterProcess(thePhotoElectricEffect, particle);
289 
290  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
291  G4LivermoreComptonModel* theLivermoreComptonModel =
293  theLivermoreComptonModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
294  theComptonScattering->AddEmModel(0, theLivermoreComptonModel);
295  ph->RegisterProcess(theComptonScattering, particle);
296 
297  G4GammaConversion* theGammaConversion = new G4GammaConversion();
298  G4LivermoreGammaConversionModel* theLivermoreGammaConversionModel =
300  theLivermoreGammaConversionModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
301  theGammaConversion->AddEmModel(0, theLivermoreGammaConversionModel);
302  ph->RegisterProcess(theGammaConversion, particle);
303 
304  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
305  G4LivermoreRayleighModel* theRayleighModel = new G4LivermoreRayleighModel();
306  theRayleighModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
307  theRayleigh->AddEmModel(0, theRayleighModel);
308  ph->RegisterProcess(theRayleigh, particle);
309  }
310 
311  // Warning : end of particles and processes are needed by EM Physics builders
312 
313  }
314 
315  // Deexcitation
316  //
319 }
320 
321 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
static G4GenericIon * GenericIonDefinition()
Definition: G4GenericIon.cc:88
#define G4DNABornIonisationModel
static G4LossTableManager * Instance()
void SetDeexcitationIgnoreCut(G4bool val)
G4VEmModel * EmModel(G4int index=1) const
void SetAuger(G4bool val)
int G4int
Definition: G4Types.hh:78
G4EmDNAPhysics_option2(G4int ver=1, const G4String &name="")
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)
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
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:731
static constexpr double keV
Definition: G4SIunits.hh:216
void SetAtomDeexcitation(G4VAtomDeexcitation *)
void SetFluo(G4bool val)
void SetStepLimitType(G4MscStepLimitType val)
G4ParticleDefinition * GetIon(const G4String &name)
#define G4_DECLARE_PHYSCONSTR_FACTORY(physics_constructor)