Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EmDNAPhysics_option5.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"
42
43#include "G4DNAIonisation.hh"
45
46#include "G4DNAExcitation.hh"
48
49#include "G4DNAAttachment.hh"
50#include "G4DNAVibExcitation.hh"
51
54
55// particles
56
57#include "G4Electron.hh"
58#include "G4Proton.hh"
59#include "G4GenericIon.hh"
60
61// Warning : the following is needed in order to use EM Physics builders
62// e+
63#include "G4Positron.hh"
65#include "G4eIonisation.hh"
66#include "G4eBremsstrahlung.hh"
68// gamma
69#include "G4Gamma.hh"
74#include "G4GammaConversion.hh"
78
79#include "G4EmParameters.hh"
80// end of warning
81
82#include "G4LossTableManager.hh"
85#include "G4BuilderType.hh"
86
87// factory
89//
91
92//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
93
95 G4VPhysicsConstructor("G4EmDNAPhysics_option5"), verbose(ver)
96{
98 param->SetDefaults();
99 param->SetFluo(true);
100 param->SetAuger(true);
101 param->SetAugerCascade(true);
102 param->SetDeexcitationIgnoreCut(true);
103 param->ActivateDNA();
104
106}
107
108//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
109
111{
112}
113
114//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
115
117{
118// bosons
120
121// leptons
124
125// baryons
127
129
130 G4DNAGenericIonsManager * genericIonsManager;
131 genericIonsManager = G4DNAGenericIonsManager::Instance();
132 genericIonsManager->GetIon("alpha++");
133 genericIonsManager->GetIon("alpha+");
134 genericIonsManager->GetIon("helium");
135 genericIonsManager->GetIon("hydrogen");
136
137}
138
139//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
140
142{
143 if(verbose > 1)
144 {
145 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
146 }
148
149 auto myParticleIterator=GetParticleIterator();
150 myParticleIterator->reset();
151 while( (*myParticleIterator)() )
152 {
153 G4ParticleDefinition* particle = myParticleIterator->value();
154 G4String particleName = particle->GetParticleName();
155
156 if (particleName == "e-")
157 {
158 // *** Solvation ***
159
160 G4DNAElectronSolvation* solvation =
161 new G4DNAElectronSolvation("e-_G4DNAElectronSolvation");
163 therm->SetHighEnergyLimit(10.*eV); // limit of the Uehara's model
164 solvation->SetEmModel(therm);
165 ph->RegisterProcess(solvation, particle);
166
167 // *** Elastic scattering ***
168 G4DNAElastic* theDNAElasticProcess = new G4DNAElastic("e-_G4DNAElastic");
170 mod->SelectFasterComputation(true);
171 theDNAElasticProcess->SetEmModel(mod);
172 ph->RegisterProcess(theDNAElasticProcess, particle);
173
174 // *** Excitation ***
175 G4DNAExcitation* theDNAExcitationProcess = new G4DNAExcitation("e-_G4DNAExcitation");
176 theDNAExcitationProcess->SetEmModel(new G4DNAEmfietzoglouExcitationModel());
177 ph->RegisterProcess(theDNAExcitationProcess, particle);
178
179 // *** Ionisation ***
180 G4DNAIonisation* theDNAIonisationProcess = new G4DNAIonisation("e-_G4DNAIonisation");
182 theDNAIonisationProcess->SetEmModel(modE);
183 modE->SelectFasterComputation(true);
184 ph->RegisterProcess(theDNAIonisationProcess, particle);
185
186 // *** Vibrational excitation ***
187 //ph->RegisterProcess(new G4DNAVibExcitation("e-_G4DNAVibExcitation"), particle);
188
189 // *** Attachment ***
190 //ph->RegisterProcess(new G4DNAAttachment("e-_G4DNAAttachment"), particle);
191
192 } else if ( particleName == "proton" ) {
193
194 ph->RegisterProcess(new G4DNAElastic("proton_G4DNAElastic"), particle);
195
196 ph->RegisterProcess(new G4DNAExcitation("proton_G4DNAExcitation"), particle);
197
198 G4DNAIonisation* theDNAIonisationProcess = new G4DNAIonisation("proton_G4DNAIonisation");
199
201 mod1->SetLowEnergyLimit(0*eV);
202 mod1->SetHighEnergyLimit(500*keV);
203
205 mod2->SetLowEnergyLimit(500*keV);
206 mod2->SetHighEnergyLimit(100*MeV);
207 mod2->SelectFasterComputation(true);
208
209 theDNAIonisationProcess->SetEmModel(mod1);
210 theDNAIonisationProcess->SetEmModel(mod2);
211
212 ph->RegisterProcess(theDNAIonisationProcess, particle);
213
214 ph->RegisterProcess(new G4DNAChargeDecrease("proton_G4DNAChargeDecrease"), particle);
215
216 } else if ( particleName == "hydrogen" ) {
217
218 ph->RegisterProcess(new G4DNAElastic("hydrogen_G4DNAElastic"), particle);
219
220 ph->RegisterProcess(new G4DNAExcitation("hydrogen_G4DNAExcitation"), particle);
221
222 G4DNAIonisation* theDNAIonisationProcess = new G4DNAIonisation("hydrogen_G4DNAIonisation");
223 theDNAIonisationProcess->SetEmModel(new G4DNARuddIonisationExtendedModel());
224 ph->RegisterProcess(theDNAIonisationProcess, particle);
225
226 ph->RegisterProcess(new G4DNAChargeIncrease("hydrogen_G4DNAChargeIncrease"), particle);
227
228 } else if ( particleName == "alpha" ) {
229
230 ph->RegisterProcess(new G4DNAElastic("alpha_G4DNAElastic"), particle);
231
232 ph->RegisterProcess(new G4DNAExcitation("alpha_G4DNAExcitation"), particle);
233
234 G4DNAIonisation* theDNAIonisationProcess = new G4DNAIonisation("alpha_G4DNAIonisation");
235 theDNAIonisationProcess->SetEmModel(new G4DNARuddIonisationExtendedModel());
236 ph->RegisterProcess(theDNAIonisationProcess, particle);
237
238 ph->RegisterProcess(new G4DNAChargeDecrease("alpha_G4DNAChargeDecrease"), particle);
239
240 } else if ( particleName == "alpha+" ) {
241
242 ph->RegisterProcess(new G4DNAElastic("alpha+_G4DNAElastic"), particle);
243
244 ph->RegisterProcess(new G4DNAExcitation("alpha+_G4DNAExcitation"), particle);
245
246 G4DNAIonisation* theDNAIonisationProcess = new G4DNAIonisation("alpha+_G4DNAIonisation");
247 theDNAIonisationProcess->SetEmModel(new G4DNARuddIonisationExtendedModel());
248 ph->RegisterProcess(theDNAIonisationProcess, particle);
249
250 ph->RegisterProcess(new G4DNAChargeDecrease("alpha+_G4DNAChargeDecrease"), particle);
251 ph->RegisterProcess(new G4DNAChargeIncrease("alpha+_G4DNAChargeIncrease"), particle);
252
253 } else if ( particleName == "helium" ) {
254
255 ph->RegisterProcess(new G4DNAElastic("helium_G4DNAElastic"), particle);
256
257 ph->RegisterProcess(new G4DNAExcitation("helium_G4DNAExcitation"), particle);
258
259 G4DNAIonisation* theDNAIonisationProcess = new G4DNAIonisation("helium_G4DNAIonisation");
260 theDNAIonisationProcess->SetEmModel(new G4DNARuddIonisationExtendedModel());
261 ph->RegisterProcess(theDNAIonisationProcess, particle);
262
263 ph->RegisterProcess(new G4DNAChargeIncrease("helium_G4DNAChargeIncrease"), particle);
264
265 // Extension to HZE proposed by Z. Francis
266
267 } else if ( particleName == "GenericIon" ) {
268 ph->RegisterProcess(new G4DNAIonisation("GenericIon_G4DNAIonisation"), particle);
269 }
270
271 // Warning : the following particles and processes are needed by EM Physics builders
272 // They are taken from the default Livermore Physics list
273 // These particles are currently not handled by Geant4-DNA
274
275 // e+
276
277 else if (particleName == "e+") {
278
281 G4eIonisation* eIoni = new G4eIonisation();
282 eIoni->SetStepFunction(0.2, 100*um);
283
284 ph->RegisterProcess(msc, particle);
285 ph->RegisterProcess(eIoni, particle);
286 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
287 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
288
289 } else if (particleName == "gamma") {
290
291 // photoelectric effect - Livermore model only
292 G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
293 thePhotoElectricEffect->SetEmModel(new G4LivermorePhotoElectricModel());
294 ph->RegisterProcess(thePhotoElectricEffect, particle);
295
296 // Compton scattering - Livermore model only
297 G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
298 theComptonScattering->SetEmModel(new G4LivermoreComptonModel());
299 ph->RegisterProcess(theComptonScattering, particle);
300
301 // gamma conversion - Livermore model below 80 GeV
302 G4GammaConversion* theGammaConversion = new G4GammaConversion();
303 theGammaConversion->SetEmModel(new G4LivermoreGammaConversionModel());
304 ph->RegisterProcess(theGammaConversion, particle);
305
306 // default Rayleigh scattering is Livermore
307 G4RayleighScattering* theRayleigh = new G4RayleighScattering();
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......
@ bElectromagnetic
#define G4DNABornIonisationModel
@ 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_option5(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 SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:764
void SetEmModel(G4VEmModel *, G4int index=0)
void SetStepFunction(G4double v1, G4double v2)
void SetStepLimitType(G4MscStepLimitType val)
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
const G4String & GetPhysicsName() const