Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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->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
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
130}
131
132//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
133
135{
136 if(verbose > 1) {
137 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
138 }
140
141 auto myParticleIterator=GetParticleIterator();
142 myParticleIterator->reset();
143 while( (*myParticleIterator)() )
144 {
145 G4ParticleDefinition* particle = myParticleIterator->value();
146 G4String particleName = particle->GetParticleName();
147
148 if (particleName == "e-") {
149
150 G4DNAElectronSolvation* solvation =
151 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 // photoelectric effect - Livermore model only
238 G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
239 thePhotoElectricEffect->SetEmModel(new G4LivermorePhotoElectricModel());
240 ph->RegisterProcess(thePhotoElectricEffect, particle);
241
242 // Compton scattering - Livermore model only
243 G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
244 theComptonScattering->SetEmModel(new G4LivermoreComptonModel());
245 ph->RegisterProcess(theComptonScattering, particle);
246
247 // gamma conversion - Livermore model below 80 GeV
248 G4GammaConversion* theGammaConversion = new G4GammaConversion();
249 theGammaConversion->SetEmModel(new G4LivermoreGammaConversionModel());
250 ph->RegisterProcess(theGammaConversion, particle);
251
252 // default Rayleigh scattering is Livermore
253 G4RayleighScattering* theRayleigh = new G4RayleighScattering();
254 ph->RegisterProcess(theRayleigh, particle);
255 }
256
257 // Warning : end of particles and processes are needed by EM Physics builders
258
259 }
260
261 // Deexcitation
262 //
265}
266
267//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ bElectromagnetic
@ fUseDistanceToBoundary
#define G4_DECLARE_PHYSCONSTR_FACTORY(physics_constructor)
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 SetAugerCascade(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