Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EmDNAPhysics_option8.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// S. Incerti ([email protected])
27//
28
30
31#include "G4SystemOfUnits.hh"
32
34
35// *** Processes and models for Geant4-DNA
36
38#include "G4DNAElastic.hh"
40
41#include "G4DNAIonisation.hh"
42
43#include "G4DNAExcitation.hh"
44
45#include "G4DNAAttachment.hh"
46#include "G4DNAVibExcitation.hh"
47
50
51// particles
52
53#include "G4Electron.hh"
54#include "G4Proton.hh"
55#include "G4GenericIon.hh"
56
57// Warning : the following is needed in order to use EM Physics builders
58// e+
59#include "G4Positron.hh"
61#include "G4eIonisation.hh"
62#include "G4eBremsstrahlung.hh"
64// gamma
65#include "G4Gamma.hh"
70#include "G4GammaConversion.hh"
74
75#include "G4EmParameters.hh"
76// end of warning
77
78#include "G4LossTableManager.hh"
81#include "G4BuilderType.hh"
82
83// factory
85//
87
88//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
89
91 : G4VPhysicsConstructor("G4EmDNAPhysics_option8"), verbose(ver)
92{
94 param->SetDefaults();
95 param->SetFluo(true);
96 param->SetAuger(true);
97 param->SetAugerCascade(true);
98 param->SetDeexcitationIgnoreCut(true);
99 param->ActivateDNA();
100
102}
103
104//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
105
107{}
108
109//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
110
112{
113// bosons
115
116// leptons
119
120// baryons
122
124
125 G4DNAGenericIonsManager * genericIonsManager;
126 genericIonsManager=G4DNAGenericIonsManager::Instance();
127 genericIonsManager->GetIon("alpha++");
128 genericIonsManager->GetIon("alpha+");
129 genericIonsManager->GetIon("helium");
130 genericIonsManager->GetIon("hydrogen");
131
132}
133
134//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
135
137{
138 if(verbose > 1) {
139 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
140 }
142
143 auto myParticleIterator=GetParticleIterator();
144 myParticleIterator->reset();
145 while( (*myParticleIterator)() )
146 {
147 G4ParticleDefinition* particle = myParticleIterator->value();
148 G4String particleName = particle->GetParticleName();
149
150 if (particleName == "e-") {
151
152 // *** Solvation ***
153
154 G4DNAElectronSolvation* solvation =
155 new G4DNAElectronSolvation("e-_G4DNAElectronSolvation");
157 therm->SetHighEnergyLimit(11.*eV); // limit of the CPA100 model
158 solvation->SetEmModel(therm);
159 ph->RegisterProcess(solvation, particle);
160
161 // *** Elastic scattering using CPA100 model ***
162
163 G4DNAElastic* theDNAElasticProcess = new G4DNAElastic("e-_G4DNAElastic");
164 G4VEmModel* model1 = new G4DNACPA100ElasticModel();
166
167 model1->SetActivationLowEnergyLimit(11*eV);
168 model1->SetActivationHighEnergyLimit(255.955*keV);
169 theDNAElasticProcess->SetEmModel(model1);
170
171 model2->SetActivationLowEnergyLimit(255.955*keV);
172 model2->SetActivationHighEnergyLimit(1*MeV);
173 theDNAElasticProcess->AddEmModel(2,model2);
174
175 ph->RegisterProcess(theDNAElasticProcess, particle);
176
177 // *** Excitation ***
178 ph->RegisterProcess(new G4DNAExcitation("e-_G4DNAExcitation"), particle);
179
180 // *** Ionisation ***
181 ph->RegisterProcess(new G4DNAIonisation("e-_G4DNAIonisation"), particle);
182
183 // *** Vibrational excitation ***
184 ph->RegisterProcess(new G4DNAVibExcitation("e-_G4DNAVibExcitation"), particle);
185
186 // *** Attachment ***
187 ph->RegisterProcess(new G4DNAAttachment("e-_G4DNAAttachment"), particle);
188
189 } else if ( particleName == "proton" ) {
190 ph->RegisterProcess(new G4DNAElastic("proton_G4DNAElastic"), particle);
191 ph->RegisterProcess(new G4DNAExcitation("proton_G4DNAExcitation"), particle);
192 ph->RegisterProcess(new G4DNAIonisation("proton_G4DNAIonisation"), particle);
193 ph->RegisterProcess(new G4DNAChargeDecrease("proton_G4DNAChargeDecrease"), particle);
194
195 } else if ( particleName == "hydrogen" ) {
196 ph->RegisterProcess(new G4DNAElastic("hydrogen_G4DNAElastic"), particle);
197 ph->RegisterProcess(new G4DNAExcitation("hydrogen_G4DNAExcitation"), particle);
198 ph->RegisterProcess(new G4DNAIonisation("hydrogen_G4DNAIonisation"), particle);
199 ph->RegisterProcess(new G4DNAChargeIncrease("hydrogen_G4DNAChargeIncrease"), particle);
200
201 } else if ( particleName == "alpha" ) {
202 ph->RegisterProcess(new G4DNAElastic("alpha_G4DNAElastic"), particle);
203 ph->RegisterProcess(new G4DNAExcitation("alpha_G4DNAExcitation"), particle);
204 ph->RegisterProcess(new G4DNAIonisation("alpha_G4DNAIonisation"), particle);
205 ph->RegisterProcess(new G4DNAChargeDecrease("alpha_G4DNAChargeDecrease"), particle);
206
207 } else if ( particleName == "alpha+" ) {
208 ph->RegisterProcess(new G4DNAElastic("alpha+_G4DNAElastic"), particle);
209 ph->RegisterProcess(new G4DNAExcitation("alpha+_G4DNAExcitation"), particle);
210 ph->RegisterProcess(new G4DNAIonisation("alpha+_G4DNAIonisation"), particle);
211 ph->RegisterProcess(new G4DNAChargeDecrease("alpha+_G4DNAChargeDecrease"), particle);
212 ph->RegisterProcess(new G4DNAChargeIncrease("alpha+_G4DNAChargeIncrease"), particle);
213
214 } else if ( particleName == "helium" ) {
215 ph->RegisterProcess(new G4DNAElastic("helium_G4DNAElastic"), particle);
216 ph->RegisterProcess(new G4DNAExcitation("helium_G4DNAExcitation"), particle);
217 ph->RegisterProcess(new G4DNAIonisation("helium_G4DNAIonisation"), particle);
218 ph->RegisterProcess(new G4DNAChargeIncrease("helium_G4DNAChargeIncrease"), particle);
219
220 // Extension to HZE proposed by Z. Francis
221
222 } else if ( particleName == "GenericIon" ) {
223 ph->RegisterProcess(new G4DNAIonisation("GenericIon_G4DNAIonisation"), particle);
224 }
225
226 // Warning : the following particles and processes are needed by EM Physics builders
227 // They are taken from the default Livermore Physics list
228 // These particles are currently not handled by Geant4-DNA
229
230 // e+
231
232 else if (particleName == "e+") {
233
234 // Identical to G4EmStandardPhysics_option3
235
238 G4eIonisation* eIoni = new G4eIonisation();
239 eIoni->SetStepFunction(0.2, 100*um);
240
241 ph->RegisterProcess(msc, particle);
242 ph->RegisterProcess(eIoni, particle);
243 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
244 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
245
246 } else if (particleName == "gamma") {
247
248 // photoelectric effect - Livermore model only
249 G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
250 thePhotoElectricEffect->SetEmModel(new G4LivermorePhotoElectricModel());
251 ph->RegisterProcess(thePhotoElectricEffect, particle);
252
253 // Compton scattering - Livermore model only
254 G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
255 theComptonScattering->SetEmModel(new G4LivermoreComptonModel());
256 ph->RegisterProcess(theComptonScattering, particle);
257
258 // gamma conversion - Livermore model below 80 GeV
259 G4GammaConversion* theGammaConversion = new G4GammaConversion();
260 theGammaConversion->SetEmModel(new G4LivermoreGammaConversionModel());
261 ph->RegisterProcess(theGammaConversion, particle);
262
263 // default Rayleigh scattering is Livermore
264 G4RayleighScattering* theRayleigh = new G4RayleighScattering();
265 ph->RegisterProcess(theRayleigh, particle);
266 }
267
268 // Warning : end of particles and processes are needed by EM Physics builders
269
270 }
271
272 // Deexcitation
273 //
276}
277
278//....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_option8(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 SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:778
void SetActivationHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:771
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetStepFunction(G4double v1, G4double v2)
void SetStepLimitType(G4MscStepLimitType val)
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
const G4String & GetPhysicsName() const