ECCE @ EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4EmDNAPhysics.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4EmDNAPhysics.cc
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 // add elastic scattering processes of proton, hydrogen, helium, alpha+, alpha++
27 
28 #include "G4EmDNAPhysics.hh"
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"), verbose(ver)
90 {
92  param->SetDefaults();
93  param->SetFluo(true);
94  param->SetAuger(true);
95  param->SetAugerCascade(true);
96  param->SetDeexcitationIgnoreCut(true);
97  param->ActivateDNA();
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 pParticleIterator = GetParticleIterator();
146  pParticleIterator->reset();
147  while( (*pParticleIterator)() )
148  {
149  G4ParticleDefinition* pParticle = pParticleIterator->value();
150  G4String particleName = pParticle->GetParticleName();
151 
152  if (particleName == "e-") {
153 
154  G4DNAElectronSolvation* pSolvation = new G4DNAElectronSolvation("e-_G4DNAElectronSolvation");
155  auto pSolvationModel = G4DNASolvationModelFactory::GetMacroDefinedModel();
156  pSolvationModel->SetHighEnergyLimit(7.4*eV); // limit of the Champion's model
157  pSolvation->SetEmModel(pSolvationModel);
158  pPhysicsHelper->RegisterProcess(pSolvation, pParticle);
159 
160  // *** Elastic scattering (two alternative models available) ***
161 
162  G4DNAElastic* pElasticProcess = new G4DNAElastic("e-_G4DNAElastic");
163  pElasticProcess->SetEmModel(new G4DNAChampionElasticModel());
164 
165  // or alternative model
166  //theDNAElasticProcess->SetEmModel(new G4DNAScreenedRutherfordElasticModel());
167 
168  pPhysicsHelper->RegisterProcess(pElasticProcess, pParticle);
169 
170  // *** Excitation ***
171  pPhysicsHelper->RegisterProcess(new G4DNAExcitation("e-_G4DNAExcitation"), pParticle);
172 
173  // *** Ionisation ***
174  pPhysicsHelper->RegisterProcess(new G4DNAIonisation("e-_G4DNAIonisation"), pParticle);
175 
176  // *** Vibrational excitation ***
177  pPhysicsHelper->RegisterProcess(new G4DNAVibExcitation("e-_G4DNAVibExcitation"), pParticle);
178 
179  // *** Attachment ***
180  pPhysicsHelper->RegisterProcess(new G4DNAAttachment("e-_G4DNAAttachment"), pParticle);
181 
182  } else if ( particleName == "proton" ) {
183  pPhysicsHelper->RegisterProcess(new G4DNAElastic("proton_G4DNAElastic"), pParticle);
184  pPhysicsHelper->RegisterProcess(new G4DNAExcitation("proton_G4DNAExcitation"), pParticle);
185  pPhysicsHelper->RegisterProcess(new G4DNAIonisation("proton_G4DNAIonisation"), pParticle);
186  pPhysicsHelper->RegisterProcess(new G4DNAChargeDecrease("proton_G4DNAChargeDecrease"), pParticle);
187 
188  } else if ( particleName == "hydrogen" ) {
189  pPhysicsHelper->RegisterProcess(new G4DNAElastic("hydrogen_G4DNAElastic"), pParticle);
190  pPhysicsHelper->RegisterProcess(new G4DNAExcitation("hydrogen_G4DNAExcitation"), pParticle);
191  pPhysicsHelper->RegisterProcess(new G4DNAIonisation("hydrogen_G4DNAIonisation"), pParticle);
192  pPhysicsHelper->RegisterProcess(new G4DNAChargeIncrease("hydrogen_G4DNAChargeIncrease"), pParticle);
193 
194  } else if ( particleName == "alpha" ) {
195  pPhysicsHelper->RegisterProcess(new G4DNAElastic("alpha_G4DNAElastic"), pParticle);
196  pPhysicsHelper->RegisterProcess(new G4DNAExcitation("alpha_G4DNAExcitation"), pParticle);
197  pPhysicsHelper->RegisterProcess(new G4DNAIonisation("alpha_G4DNAIonisation"), pParticle);
198  pPhysicsHelper->RegisterProcess(new G4DNAChargeDecrease("alpha_G4DNAChargeDecrease"), pParticle);
199 
200  } else if ( particleName == "alpha+" ) {
201  pPhysicsHelper->RegisterProcess(new G4DNAElastic("alpha+_G4DNAElastic"), pParticle);
202  pPhysicsHelper->RegisterProcess(new G4DNAExcitation("alpha+_G4DNAExcitation"), pParticle);
203  pPhysicsHelper->RegisterProcess(new G4DNAIonisation("alpha+_G4DNAIonisation"), pParticle);
204  pPhysicsHelper->RegisterProcess(new G4DNAChargeDecrease("alpha+_G4DNAChargeDecrease"), pParticle);
205  pPhysicsHelper->RegisterProcess(new G4DNAChargeIncrease("alpha+_G4DNAChargeIncrease"), pParticle);
206 
207  } else if ( particleName == "helium" ) {
208  pPhysicsHelper->RegisterProcess(new G4DNAElastic("helium_G4DNAElastic"), pParticle);
209  pPhysicsHelper->RegisterProcess(new G4DNAExcitation("helium_G4DNAExcitation"), pParticle);
210  pPhysicsHelper->RegisterProcess(new G4DNAIonisation("helium_G4DNAIonisation"), pParticle);
211  pPhysicsHelper->RegisterProcess(new G4DNAChargeIncrease("helium_G4DNAChargeIncrease"), pParticle);
212 
213  } else if ( particleName == "GenericIon" ) {
214  pPhysicsHelper->RegisterProcess(new G4DNAIonisation("GenericIon_G4DNAIonisation"), pParticle);
215 
216  /*
217  } else if ( particleName == "carbon" ) {
218  pPhysicsHelper->RegisterProcess(new G4DNAIonisation("carbon_G4DNAIonisation"), particle);
219 
220  } else if ( particleName == "nitrogen" ) {
221  pPhysicsHelper->RegisterProcess(new G4DNAIonisation("nitrogen_G4DNAIonisation"), particle);
222 
223  } else if ( particleName == "oxygen" ) {
224  pPhysicsHelper->RegisterProcess(new G4DNAIonisation("oxygen_G4DNAIonisation"), particle);
225 
226  } else if ( particleName == "iron" ) {
227  pPhysicsHelper->RegisterProcess(new G4DNAIonisation("iron_G4DNAIonisation"), particle);
228  */
229 
230  }
231 
232  // Warning : the following particles and processes are needed by EM Physics builders
233  // They are taken from the default Livermore Physics list
234  // These particles are currently not handled by Geant4-DNA
235 
236  // e+
237 
238  else if (particleName == "e+") {
239 
242  G4eIonisation* eIoni = new G4eIonisation();
243  eIoni->SetStepFunction(0.2, 100*um);
244 
245  pPhysicsHelper->RegisterProcess(msc, pParticle);
246  pPhysicsHelper->RegisterProcess(eIoni, pParticle);
247  pPhysicsHelper->RegisterProcess(new G4eBremsstrahlung(), pParticle);
248  pPhysicsHelper->RegisterProcess(new G4eplusAnnihilation(), pParticle);
249 
250  } else if (particleName == "gamma") {
251 
252  // photoelectric effect - Livermore model only
253  G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
254  thePhotoElectricEffect->SetEmModel(new G4LivermorePhotoElectricModel());
255  pPhysicsHelper->RegisterProcess(thePhotoElectricEffect, pParticle);
256 
257  // Compton scattering - Livermore model only
258  G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
259  theComptonScattering->SetEmModel(new G4LivermoreComptonModel());
260  pPhysicsHelper->RegisterProcess(theComptonScattering, pParticle);
261 
262  // gamma conversion - Livermore model below 80 GeV
263  G4GammaConversion* theGammaConversion = new G4GammaConversion();
264  theGammaConversion->SetEmModel(new G4LivermoreGammaConversionModel());
265  pPhysicsHelper->RegisterProcess(theGammaConversion, pParticle);
266 
267  // default Rayleigh scattering is Livermore
268  G4RayleighScattering* theRayleigh = new G4RayleighScattering();
269  pPhysicsHelper->RegisterProcess(theRayleigh, pParticle);
270  }
271  // Warning : end of particles and processes are needed by EM Physics builders
272  }
273 
274  // Deexcitation
275  //
278 }
279 
280 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......