Geant4-11
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"
41
42#include "G4DNAExcitation.hh"
43#include "G4DNAAttachment.hh"
44#include "G4DNAVibExcitation.hh"
45#include "G4DNAIonisation.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"
59#include "G4eIonisation.hh"
60#include "G4eBremsstrahlung.hh"
62// gamma
63#include "G4Gamma.hh"
68#include "G4GammaConversion.hh"
72
73#include "G4EmParameters.hh"
74// end of warning
75
76#include "G4LossTableManager.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->SetDeexcitationIgnoreCut(true);
96 param->ActivateDNA();
97
99}
100
101//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
102
104{}
105
106//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
107
109{
110// bosons
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");
152 therm->SetHighEnergyLimit(7.4*eV); // limit of the Champion's model
153 solvation->SetEmModel(therm);
154 ph->RegisterProcess(solvation, particle);
155
156 // *** Elastic scattering (two alternative models available) ***
157
158 G4DNAElastic* theDNAElasticProcess = new G4DNAElastic("e-_G4DNAElastic");
159 theDNAElasticProcess->SetEmModel(new G4DNAChampionElasticModel());
160
161 // or alternative model
162 //theDNAElasticProcess->SetEmModel(new G4DNAScreenedRutherfordElasticModel());
163
164 ph->RegisterProcess(theDNAElasticProcess, particle);
165
166 // *** Excitation ***
167 ph->RegisterProcess(new G4DNAExcitation("e-_G4DNAExcitation"), particle);
168
169 // *** Ionisation ***
170 ph->RegisterProcess(new G4DNAIonisation("e-_G4DNAIonisation"), particle);
171
172 // *** Vibrational excitation ***
173 ph->RegisterProcess(new G4DNAVibExcitation("e-_G4DNAVibExcitation"), particle);
174
175 // *** Attachment ***
176 ph->RegisterProcess(new G4DNAAttachment("e-_G4DNAAttachment"), particle);
177
178 } else if ( particleName == "proton" ) {
179 ph->RegisterProcess(new G4DNAElastic("proton_G4DNAElastic"), particle);
180 ph->RegisterProcess(new G4DNAExcitation("proton_G4DNAExcitation"), particle);
181 ph->RegisterProcess(new G4DNAIonisation("proton_G4DNAIonisation"), particle);
182 ph->RegisterProcess(new G4DNAChargeDecrease("proton_G4DNAChargeDecrease"), particle);
183
184 } else if ( particleName == "hydrogen" ) {
185 ph->RegisterProcess(new G4DNAElastic("hydrogen_G4DNAElastic"), particle);
186 ph->RegisterProcess(new G4DNAExcitation("hydrogen_G4DNAExcitation"), particle);
187 ph->RegisterProcess(new G4DNAIonisation("hydrogen_G4DNAIonisation"), particle);
188 ph->RegisterProcess(new G4DNAChargeIncrease("hydrogen_G4DNAChargeIncrease"), particle);
189
190 } else if ( particleName == "alpha" ) {
191 ph->RegisterProcess(new G4DNAElastic("alpha_G4DNAElastic"), particle);
192 ph->RegisterProcess(new G4DNAExcitation("alpha_G4DNAExcitation"), particle);
193 ph->RegisterProcess(new G4DNAIonisation("alpha_G4DNAIonisation"), particle);
194 ph->RegisterProcess(new G4DNAChargeDecrease("alpha_G4DNAChargeDecrease"), 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 ph->RegisterProcess(new G4DNAChargeIncrease("alpha+_G4DNAChargeIncrease"), particle);
202
203 } else if ( particleName == "helium" ) {
204 ph->RegisterProcess(new G4DNAElastic("helium_G4DNAElastic"), particle);
205 ph->RegisterProcess(new G4DNAExcitation("helium_G4DNAExcitation"), particle);
206 ph->RegisterProcess(new G4DNAIonisation("helium_G4DNAIonisation"), particle);
207 ph->RegisterProcess(new G4DNAChargeIncrease("helium_G4DNAChargeIncrease"), particle);
208
209 } else if ( particleName == "GenericIon" ) {
210 ph->RegisterProcess(new G4DNAIonisation("GenericIon_G4DNAIonisation"), particle);
211
212 }
213
214 // Warning : the following particles and processes are needed by EM Physics builders
215 // They are taken from the default Livermore Physics list
216 // These particles are currently not handled by Geant4-DNA
217
218 // e+
219
220 else if (particleName == "e+") {
221
222 // Identical to G4EmStandardPhysics_option3
223
226 G4eIonisation* eIoni = new G4eIonisation();
227 eIoni->SetStepFunction(0.2, 100*um);
228
229 ph->RegisterProcess(msc, particle);
230 ph->RegisterProcess(eIoni, particle);
231 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
232 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
233
234 } else if (particleName == "gamma") {
235
236 // photoelectric effect - Livermore model only
237 G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
238 thePhotoElectricEffect->SetEmModel(new G4LivermorePhotoElectricModel());
239 ph->RegisterProcess(thePhotoElectricEffect, particle);
240
241 // Compton scattering - Livermore model only
242 G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
243 theComptonScattering->SetEmModel(new G4LivermoreComptonModel());
244 ph->RegisterProcess(theComptonScattering, particle);
245
246 // gamma conversion - Livermore model below 80 GeV
247 G4GammaConversion* theGammaConversion = new G4GammaConversion();
248 theGammaConversion->SetEmModel(new G4LivermoreGammaConversionModel());
249 ph->RegisterProcess(theGammaConversion, particle);
250
251 // default Rayleigh scattering is Livermore
252 G4RayleighScattering* theRayleigh = new G4RayleighScattering();
253 ph->RegisterProcess(theRayleigh, particle);
254 }
255
256 // Warning : end of particles and processes are needed by EM Physics builders
257
258 }
259
260 // Deexcitation
261 //
264}
265
266//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ bElectromagnetic
G4_DECLARE_PHYSCONSTR_FACTORY(G4EmDNAPhysics_option3)
@ fUseDistanceToBoundary
static constexpr double eV
Definition: G4SIunits.hh:201
static constexpr double um
Definition: G4SIunits.hh:93
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4DNAGenericIonsManager * Instance(void)
G4ParticleDefinition * GetIon(const G4String &name)
static G4VEmModel * GetMacroDefinedModel()
One step thermalization model can be chosen via macro using /process/dna/e-SolvationSubType Ritchie19...
static G4Electron * Electron()
Definition: G4Electron.cc:93
G4EmDNAPhysics_option3(G4int ver=1, const G4String &name="")
static G4EmParameters * Instance()
void SetDeexcitationIgnoreCut(G4bool val)
void SetFluo(G4bool val)
void SetAuger(G4bool val)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
static G4GenericIon * GenericIonDefinition()
Definition: G4GenericIon.cc:87
void SetAtomDeexcitation(G4VAtomDeexcitation *)
static G4LossTableManager * Instance()
const G4String & GetParticleName() const
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4PhysicsListHelper * GetPhysicsListHelper()
static G4Positron * Positron()
Definition: G4Positron.cc:93
static G4Proton * Proton()
Definition: G4Proton.cc:92
void SetEmModel(G4VEmModel *, G4int index=0)
void SetStepFunction(G4double v1, G4double v2)
void SetStepLimitType(G4MscStepLimitType val)
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
const G4String & GetPhysicsName() const